Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

PSF model #2643

Open
wants to merge 24 commits into
base: main
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from 22 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions docs/changes/2643.feature.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@
Adds psf models to ``ctapipe.instrument.optics`` with the parent class ``PSFModel`` and a psf model based on pure coma aberration called ``ComaModel``
- The function ``PSFModel.pdf`` gives the value of the PSF in a given location
11 changes: 11 additions & 0 deletions docs/references.bib
Original file line number Diff line number Diff line change
Expand Up @@ -137,6 +137,17 @@ @article{bright-star-catalog
url={https://heasarc.gsfc.nasa.gov/W3Browse/star-catalog/bsc5p.html},
}


@article{startracker,
author={Abe, K. and others},
title={Star tracking for pointing determination of Imaging Atmospheric Cherenkov Telescopes - Application to the Large-Sized Telescope of the Cherenkov Telescope Array},
doi={10.1051/0004-6361/202347128},
journal={A\&A},
year={2023},
volume={679},
pages={A90},
}

@manual{ctao-software-standards,
title={Software Programming Standards},
number={CTA-STD-OSO-000000-0001},
Expand Down
165 changes: 165 additions & 0 deletions src/ctapipe/instrument/optics.py
Original file line number Diff line number Diff line change
Expand Up @@ -3,11 +3,17 @@
"""

import logging
from abc import abstractmethod
from enum import Enum, auto, unique

import astropy.units as u
import numpy as np
from astropy.table import QTable
from erfa.ufunc import p2s as cartesian_to_spherical
from scipy.stats import laplace, laplace_asymmetric

from ctapipe.core import TelescopeComponent
from ctapipe.core.traits import List

from ..compat import StrEnum
from ..utils import get_table_dataset
Expand All @@ -18,6 +24,8 @@
__all__ = [
"OpticsDescription",
"FocalLengthKind",
"PSFModel",
"ComaModel",
]


Expand Down Expand Up @@ -263,3 +271,160 @@ def __repr__(self):

def __str__(self):
return self.name


class PSFModel(TelescopeComponent):
"""
Base component to describe image distortion due to the optics of the different cameras.
"""

@abstractmethod
def pdf(self, r, f, r0, f0, *args):
"""
Calculates the value of the psf at a given location

Parameters
----------
r : float
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

These docstrings are confusing to me. I don't know how to understand the "from where / at where" distinction.

One is the position of the point source and one the position where the PSF is evaluated, right?

Also: We don't really have anything that is using radial coordinates in the camera plane yet.

For the API, I think I'd prefer taking CameraFrame coordinate instances or x / y values. The conversions to polar coordinates can then be done for methods where this is needed for evaluation.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

i changed the API to use carthesian coordinates and fixed the docustrings to be more clear. I also changed the parameters to round numbers that give reasonable values for some LST-sized camera.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Here it still has r,f,r0,f0

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

And if meters are required, then this should probably use u.quantity_input(..) with astropy units

distance to the center of the camera in meters, location at where the PSF is evaluated
f : float
polar angle in radians, location at where the PSF is evaluated
r0 : float
distance to the center of the camera in meters, location from where the PSF is evaluated
f0 : float
polar angle in radians, location from where the PSF is evaluated
Returns
----------
psf : float
value of the PSF at the specified location
"""
pass


class ComaModel(PSFModel):
r"""
PSF model, describing pure coma aberrations PSF effect
asymmetry_params describe the dependency of the psf on the distance to the center of the camera
Used to calculate a pdf asymmetry parameter K of the asymmetric radial laplacian of the psf as
a function of the distance r to the optical axis

.. math:: K(r) = 1 - asym_0 \tanh(asym_1 r) - asym_2 r

radial_scale_params describes the dependency of the radial scale on the distance to the center of the camera
Used to calculate width Sr of the asymmetric radial laplacian in the PSF as a of function the distance r to the optical axis

.. math:: S_{R}(r) & = b_1 + b_2\,r + b_3\,r^2 + b_4\,r^3

az_scale_params Describes the dependency of the azimuthal scale on the distance to the center of the camera
Used to calculate the width Sf of the azimuthal laplacian in the PSF as a function of the angle :math:`phi`

.. math:: S_{\phi}(r) & = a_1\,\exp{(-a_2\,r)}+\frac{a_3}{a_3+r}

for reference see :cite:p:`startracker`
"""

asymmetry_params = List(
help=(
"Describes the dependency of the psf on the distance"
"to the center of the camera. Used to calculate a pdf"
"asymmetry parameter K of the asymmetric radial laplacian"
"of the psf as a function of the distance r to the optical axis"
)
).tag(config=True)

radial_scale_params = List(
help=(
"Describes the dependency of the radial scale on the"
"distance to the center of the camera Used to calculate"
"width Sr of the asymmetric radial laplacian in the PSF"
"as a of function the distance r to the optical axis"
)
).tag(config=True)

az_scale_params = List(
help=(
"Describes the dependency of the azimuthal scale on the"
"distance to the center of the camera. Used to calculate"
"the the width Sf of the azimuthal laplacian in the PSF"
"as a function of the angle phi"
)
).tag(config=True)

def __init__(
self,
subarray,
**kwargs,
):
r"""
PSF model, describing purely the effect of coma aberration on the PSF
uses an asymmetric laplacian for the radial part

.. math:: f_{R}(r, K) = \begin{cases}\frac{1}{S_{R}(K+K^{-1})}e^{-K\frac{r-L}{S_{R}}}, r\ge L\\ \frac{1}{S_{R}(K+K^{-1})}e^{\frac{r-L}{KS_{R}}}, r < L\end{cases}

and a symmetric laplacian in azimuthal direction

.. math:: f_{\Phi}(\phi) = \frac{1}{2S_\phi}e^{-|\frac{\phi-\phi_0}{S_\phi}|}

Parameters
----------
subarray: ctapipe.instrument.SubarrayDescription
Description of the subarray.
"""
super().__init__(subarray=subarray, **kwargs)
self.check_model_parameters()

def k_func(self, x):
return (
1
- self.asymmetry_params[0] * np.tanh(self.asymmetry_params[1] * x)
- self.asymmetry_params[2] * x
)

def sr_func(self, x):
return np.polyval(self.radial_scale_params, x)

def sf_func(self, x):
return self.az_scale_params[0] * np.exp(
-self.az_scale_params[1] * x
) + self.az_scale_params[2] / (self.az_scale_params[2] + x)

def pdf(self, x, y, x0, y0):
"""
Calculates the value of the psf at a given location

Parameters
----------
x : float
x-coordinate of the point on the focal plane where the psf is evaluated in meters
y : float
y-coordinate of the point on the focal plane where the psf is evaluated in meters
x0 : float
x-coordinate of the point source on the focal plane in meters
y0 : float
y-coordinate of the point source on the focal plane in meters
Returns
----------
psf : float
value of the PSF at the specified location with the specified position of the point source
"""
f, _, r = cartesian_to_spherical((x, y, 0.0))
f0, _, r0 = cartesian_to_spherical((x0, y0, 0.0))
k = self.k_func(r0)
sr = self.sr_func(r0)
sf = self.sf_func(r0)
self.radial_pdf_params = (k, r0, sr)
self.azimuthal_pdf_params = (f0, sf)

return laplace_asymmetric.pdf(r, *self.radial_pdf_params) * laplace.pdf(
f, *self.azimuthal_pdf_params
)

def check_model_parameters(self):
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This should use the traitlets.validate decorator, one for each traitlet instead of being called manually:

See https://traitlets.readthedocs.io/en/stable/api.html#validating-proposed-changes

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sounds good. I will get to it in about an hour.

if not (
len(self.asymmetry_params) == 3
and len(self.radial_scale_params) == 4
and len(self.az_scale_params) == 3
):
raise ValueError(
"asymmetry_params and az_scale_params needs to have length 3 and radial_scale_params length 4"
)
37 changes: 37 additions & 0 deletions src/ctapipe/instrument/tests/test_psf_model.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,37 @@
"""
This module contains the ctapipe.image.psf_model unit tests
"""
import numpy as np
import pytest

from ctapipe.instrument.optics import PSFModel


@pytest.fixture(scope="session")
def coma_psf(example_subarray):
psf = PSFModel.from_name(
"ComaModel",
subarray=example_subarray,
asymmetry_params=[0.5, 10, 0.15],
radial_scale_params=[0.015, -0.1, 0.06, 0.03],
az_scale_params=[0.25, 7.5, 0.02],
)
return psf


def test_psf(example_subarray):
with pytest.raises(
ValueError,
match="asymmetry_params and az_scale_params needs to have length 3 and radial_scale_params length 4",
):
PSFModel.from_name(
mexanick marked this conversation as resolved.
Show resolved Hide resolved
"ComaModel",
subarray=example_subarray,
asymmetry_params=[0.0, 0.0, 0.0],
radial_scale_params=[0.0, 0.0, 0.0],
az_scale_params=[0.0],
)


def test_asymptotic_behavior(coma_psf):
assert np.isclose(coma_psf.pdf(10.0, 0.0, 1.0, 0.0), 0.0)
Loading