LSST Applications g0603fd7c41+f1f8eaba91,g124d44cf3d+ce19972735,g180d380827+c1373eaf06,g1afd7665f7+eb25d4c773,g2079a07aa2+86d27d4dc4,g2305ad1205+aa3c8c93b6,g2bbee38e9b+44a02a0554,g337abbeb29+44a02a0554,g33d1c0ed96+44a02a0554,g3a166c0a6a+44a02a0554,g3d1719c13e+a4710a6d26,g487adcacf7+e387efc8c5,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+7f57e6be76,g858d7b2824+a4710a6d26,g991b906543+a4710a6d26,g99cad8db69+832a1c95fd,g9b9dfce982+e7b986f76c,g9ddcbc5298+9a081db1e4,ga1e77700b3+03d07e1c1f,gb0e22166c9+60f28cb32d,gb23b769143+a4710a6d26,gb3a676b8dc+e2510deafe,gba4ed39666+c2a2e4ac27,gbb8dafda3b+201573ceae,gbd998247f1+585e252eca,gc120e1dc64+7fb97cd961,gc28159a63d+44a02a0554,gc3e9b769f7+20d5ea8805,gcf0d15dbbd+e7b986f76c,gdaeeff99f8+f9a426f77a,ge79ae78c31+44a02a0554,ged0e8a7f67+8df1cf93fe,gee10cc3b42+585e252eca,w.2024.18
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
lsst.meas.algorithms.maskStreaks.LineProfile Class Reference

Public Member Functions

 __init__ (self, data, weights, line=None)
 
 setLineMask (self, line)
 
 makeProfile (self, line, fitFlux=True)
 
 fit (self, dChi2Tol=0.1, maxIter=100, log=None)
 

Public Attributes

 data
 
 weights
 
 mask
 
 lineMask
 
 lineMaskSize
 

Protected Member Functions

 _makeMaskedProfile (self, line, fitFlux=True)
 
 _lineChi2 (self, line, grad=True)
 

Protected Attributes

 _ymax
 
 _xmax
 
 _dtype
 
 _rhoMax
 
 _xmesh
 
 _ymesh
 
 _initLine
 
 _maskData
 
 _maskWeights
 
 _mxmesh
 
 _mymesh
 

Detailed Description

Construct and/or fit a model for a linear streak.

This assumes a simple model for a streak, in which the streak
follows a straight line in pixels space, with a Moffat-shaped profile. The
model is fit to data using a Newton-Raphson style minimization algorithm.
The initial guess for the line parameters is assumed to be fairly accurate,
so only a narrow band of pixels around the initial line estimate is used in
fitting the model, which provides a significant speed-up over using all the
data. The class can also be used just to construct a model for the data with
a line following the given coordinates.

Parameters
----------
data : `np.ndarray`
    2d array of data.
weights : `np.ndarray`
    2d array of weights.
line : `Line`, optional
    Guess for position of line. Data far from line guess is masked out.
    Defaults to None, in which case only data with `weights` = 0 is masked
    out.

Definition at line 106 of file maskStreaks.py.

Constructor & Destructor Documentation

◆ __init__()

lsst.meas.algorithms.maskStreaks.LineProfile.__init__ ( self,
data,
weights,
line = None )

Definition at line 130 of file maskStreaks.py.

130 def __init__(self, data, weights, line=None):
131 self.data = data
132 self.weights = weights
133 self._ymax, self._xmax = data.shape
134 self._dtype = data.dtype
135 xrange = np.arange(self._xmax) - self._xmax / 2.
136 yrange = np.arange(self._ymax) - self._ymax / 2.
137 self._rhoMax = ((0.5 * self._ymax)**2 + (0.5 * self._xmax)**2)**0.5
138 self._xmesh, self._ymesh = np.meshgrid(xrange, yrange)
139 self.mask = (weights != 0)
140
141 self._initLine = line
142 self.setLineMask(line)
143

Member Function Documentation

◆ _lineChi2()

lsst.meas.algorithms.maskStreaks.LineProfile._lineChi2 ( self,
line,
grad = True )
protected
Construct the chi2 between the data and the model.

Parameters
----------
line : `Line`
    `Line` parameters for which to build model and calculate chi2.
grad : `bool`, optional
    Whether or not to return the gradient and hessian.

Returns
-------
reducedChi : `float`
    Reduced chi2 of the model.
reducedDChi : `np.ndarray`
    Derivative of the chi2 with respect to rho, theta, invSigma.
reducedHessianChi : `np.ndarray`
    Hessian of the chi2 with respect to rho, theta, invSigma.

Definition at line 245 of file maskStreaks.py.

245 def _lineChi2(self, line, grad=True):
246 """Construct the chi2 between the data and the model.
247
248 Parameters
249 ----------
250 line : `Line`
251 `Line` parameters for which to build model and calculate chi2.
252 grad : `bool`, optional
253 Whether or not to return the gradient and hessian.
254
255 Returns
256 -------
257 reducedChi : `float`
258 Reduced chi2 of the model.
259 reducedDChi : `np.ndarray`
260 Derivative of the chi2 with respect to rho, theta, invSigma.
261 reducedHessianChi : `np.ndarray`
262 Hessian of the chi2 with respect to rho, theta, invSigma.
263 """
264 # Calculate chi2
265 model, dModel = self._makeMaskedProfile(line)
266 chi2 = (self._maskWeights * (self._maskData - model)**2).sum()
267 if not grad:
268 return chi2.sum() / self.lineMaskSize
269
270 # Calculate derivative and Hessian of chi2
271 derivChi2 = ((-2 * self._maskWeights * (self._maskData - model))[None, :] * dModel).sum(axis=1)
272 hessianChi2 = (2 * self._maskWeights * dModel[:, None, :] * dModel[None, :, :]).sum(axis=2)
273
274 reducedChi = chi2 / self.lineMaskSize
275 reducedDChi = derivChi2 / self.lineMaskSize
276 reducedHessianChi = hessianChi2 / self.lineMaskSize
277 return reducedChi, reducedDChi, reducedHessianChi
278

◆ _makeMaskedProfile()

lsst.meas.algorithms.maskStreaks.LineProfile._makeMaskedProfile ( self,
line,
fitFlux = True )
protected
Construct the line model in the masked region and calculate its
derivatives.

Parameters
----------
line : `Line`
    Parameters of line profile for which to make profile in the masked
    region.
fitFlux : `bool`
    Fit the amplitude of the line profile to the data.

Returns
-------
model : `np.ndarray`
    Model in the masked region.
dModel : `np.ndarray`
    Derivative of the model in the masked region.

Definition at line 167 of file maskStreaks.py.

167 def _makeMaskedProfile(self, line, fitFlux=True):
168 """Construct the line model in the masked region and calculate its
169 derivatives.
170
171 Parameters
172 ----------
173 line : `Line`
174 Parameters of line profile for which to make profile in the masked
175 region.
176 fitFlux : `bool`
177 Fit the amplitude of the line profile to the data.
178
179 Returns
180 -------
181 model : `np.ndarray`
182 Model in the masked region.
183 dModel : `np.ndarray`
184 Derivative of the model in the masked region.
185 """
186 invSigma = line.sigma**-1
187 # Calculate distance between pixels and line
188 radtheta = np.deg2rad(line.theta)
189 costheta = np.cos(radtheta)
190 sintheta = np.sin(radtheta)
191 distance = (costheta * self._mxmesh + sintheta * self._mymesh - line.rho)
192 distanceSquared = distance**2
193
194 # Calculate partial derivatives of distance
195 drad = np.pi / 180
196 dDistanceSqdRho = 2 * distance * (-np.ones_like(self._mxmesh))
197 dDistanceSqdTheta = (2 * distance * (-sintheta * self._mxmesh + costheta * self._mymesh) * drad)
198
199 # Use pixel-line distances to make Moffat profile
200 profile = (1 + distanceSquared * invSigma**2)**-2.5
201 dProfile = -2.5 * (1 + distanceSquared * invSigma**2)**-3.5
202
203 if fitFlux:
204 # Calculate line flux from profile and data
205 flux = ((self._maskWeights * self._maskData * profile).sum()
206 / (self._maskWeights * profile**2).sum())
207 else:
208 # Approximately normalize the line
209 flux = invSigma**-1
210 if np.isnan(flux):
211 flux = 0
212
213 model = flux * profile
214
215 # Calculate model derivatives
216 fluxdProfile = flux * dProfile
217 fluxdProfileInvSigma = fluxdProfile * invSigma**2
218 dModeldRho = fluxdProfileInvSigma * dDistanceSqdRho
219 dModeldTheta = fluxdProfileInvSigma * dDistanceSqdTheta
220 dModeldInvSigma = fluxdProfile * distanceSquared * 2 * invSigma
221
222 dModel = np.array([dModeldRho, dModeldTheta, dModeldInvSigma])
223 return model, dModel
224

◆ fit()

lsst.meas.algorithms.maskStreaks.LineProfile.fit ( self,
dChi2Tol = 0.1,
maxIter = 100,
log = None )
Perform Newton-Raphson minimization to find line parameters.

This method takes advantage of having known derivative and Hessian of
the multivariate function to quickly and efficiently find the minimum.
This is more efficient than the scipy implementation of the Newton-
Raphson method, which doesn't take advantage of the Hessian matrix. The
method here also performs a line search in the direction of the steepest
derivative at each iteration, which reduces the number of iterations
needed.

Parameters
----------
dChi2Tol : `float`, optional
    Change in Chi2 tolerated for fit convergence.
maxIter : `int`, optional
    Maximum number of fit iterations allowed. The fit should converge in
    ~10 iterations, depending on the value of dChi2Tol, but this
    maximum provides a backup.
log : `lsst.utils.logging.LsstLogAdapter`, optional
    Logger to use for reporting more details for fitting failures.

Returns
-------
outline : `np.ndarray`
    Coordinates and inverse width of fit line.
chi2 : `float`
    Reduced Chi2 of model fit to data.
fitFailure : `bool`
    Boolean where `False` corresponds to a successful  fit.

Definition at line 279 of file maskStreaks.py.

279 def fit(self, dChi2Tol=0.1, maxIter=100, log=None):
280 """Perform Newton-Raphson minimization to find line parameters.
281
282 This method takes advantage of having known derivative and Hessian of
283 the multivariate function to quickly and efficiently find the minimum.
284 This is more efficient than the scipy implementation of the Newton-
285 Raphson method, which doesn't take advantage of the Hessian matrix. The
286 method here also performs a line search in the direction of the steepest
287 derivative at each iteration, which reduces the number of iterations
288 needed.
289
290 Parameters
291 ----------
292 dChi2Tol : `float`, optional
293 Change in Chi2 tolerated for fit convergence.
294 maxIter : `int`, optional
295 Maximum number of fit iterations allowed. The fit should converge in
296 ~10 iterations, depending on the value of dChi2Tol, but this
297 maximum provides a backup.
298 log : `lsst.utils.logging.LsstLogAdapter`, optional
299 Logger to use for reporting more details for fitting failures.
300
301 Returns
302 -------
303 outline : `np.ndarray`
304 Coordinates and inverse width of fit line.
305 chi2 : `float`
306 Reduced Chi2 of model fit to data.
307 fitFailure : `bool`
308 Boolean where `False` corresponds to a successful fit.
309 """
310 # Do minimization on inverse of sigma to simplify derivatives:
311 x = np.array([self._initLine.rho, self._initLine.theta, self._initLine.sigma**-1])
312
313 dChi2 = 1
314 iter = 0
315 oldChi2 = 0
316 fitFailure = False
317
318 def line_search(c, dx):
319 testx = x - c * dx
320 testLine = Line(testx[0], testx[1], testx[2]**-1)
321 return self._lineChi2(testLine, grad=False)
322
323 while abs(dChi2) > dChi2Tol:
324 line = Line(x[0], x[1], x[2]**-1)
325 chi2, b, A = self._lineChi2(line)
326 if chi2 == 0:
327 break
328 if not np.isfinite(A).all():
329 fitFailure = True
330 if log is not None:
331 log.warning("Hessian matrix has non-finite elements.")
332 break
333 dChi2 = oldChi2 - chi2
334 try:
335 cholesky = scipy.linalg.cho_factor(A)
336 except np.linalg.LinAlgError:
337 fitFailure = True
338 if log is not None:
339 log.warning("Hessian matrix is not invertible.")
340 break
341 dx = scipy.linalg.cho_solve(cholesky, b)
342
343 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05)
344 x -= factor * dx
345 if (abs(x[0]) > 1.5 * self._rhoMax) or (iter > maxIter):
346 fitFailure = True
347 break
348 oldChi2 = chi2
349 iter += 1
350
351 outline = Line(x[0], x[1], abs(x[2])**-1)
352
353 return outline, chi2, fitFailure
354
355

◆ makeProfile()

lsst.meas.algorithms.maskStreaks.LineProfile.makeProfile ( self,
line,
fitFlux = True )
Construct the line profile model.

Parameters
----------
line : `Line`
    Parameters of the line profile to model.
fitFlux : `bool`, optional
    Fit the amplitude of the line profile to the data.

Returns
-------
finalModel : `np.ndarray`
    Model for line profile.

Definition at line 225 of file maskStreaks.py.

225 def makeProfile(self, line, fitFlux=True):
226 """Construct the line profile model.
227
228 Parameters
229 ----------
230 line : `Line`
231 Parameters of the line profile to model.
232 fitFlux : `bool`, optional
233 Fit the amplitude of the line profile to the data.
234
235 Returns
236 -------
237 finalModel : `np.ndarray`
238 Model for line profile.
239 """
240 model, _ = self._makeMaskedProfile(line, fitFlux=fitFlux)
241 finalModel = np.zeros((self._ymax, self._xmax), dtype=self._dtype)
242 finalModel[self.lineMask] = model
243 return finalModel
244

◆ setLineMask()

lsst.meas.algorithms.maskStreaks.LineProfile.setLineMask ( self,
line )
Set mask around the image region near the line.

Parameters
----------
line : `Line`
    Parameters of line in the image.

Definition at line 144 of file maskStreaks.py.

144 def setLineMask(self, line):
145 """Set mask around the image region near the line.
146
147 Parameters
148 ----------
149 line : `Line`
150 Parameters of line in the image.
151 """
152 if line:
153 # Only fit pixels within 5 sigma of the estimated line
154 radtheta = np.deg2rad(line.theta)
155 distance = (np.cos(radtheta) * self._xmesh + np.sin(radtheta) * self._ymesh - line.rho)
156 m = (abs(distance) < 5 * line.sigma)
157 self.lineMask = self.mask & m
158 else:
159 self.lineMask = np.copy(self.mask)
160
161 self.lineMaskSize = self.lineMask.sum()
162 self._maskData = self.data[self.lineMask]
163 self._maskWeights = self.weights[self.lineMask]
164 self._mxmesh = self._xmesh[self.lineMask]
165 self._mymesh = self._ymesh[self.lineMask]
166

Member Data Documentation

◆ _dtype

lsst.meas.algorithms.maskStreaks.LineProfile._dtype
protected

Definition at line 134 of file maskStreaks.py.

◆ _initLine

lsst.meas.algorithms.maskStreaks.LineProfile._initLine
protected

Definition at line 141 of file maskStreaks.py.

◆ _maskData

lsst.meas.algorithms.maskStreaks.LineProfile._maskData
protected

Definition at line 162 of file maskStreaks.py.

◆ _maskWeights

lsst.meas.algorithms.maskStreaks.LineProfile._maskWeights
protected

Definition at line 163 of file maskStreaks.py.

◆ _mxmesh

lsst.meas.algorithms.maskStreaks.LineProfile._mxmesh
protected

Definition at line 164 of file maskStreaks.py.

◆ _mymesh

lsst.meas.algorithms.maskStreaks.LineProfile._mymesh
protected

Definition at line 165 of file maskStreaks.py.

◆ _rhoMax

lsst.meas.algorithms.maskStreaks.LineProfile._rhoMax
protected

Definition at line 137 of file maskStreaks.py.

◆ _xmax

lsst.meas.algorithms.maskStreaks.LineProfile._xmax
protected

Definition at line 133 of file maskStreaks.py.

◆ _xmesh

lsst.meas.algorithms.maskStreaks.LineProfile._xmesh
protected

Definition at line 138 of file maskStreaks.py.

◆ _ymax

lsst.meas.algorithms.maskStreaks.LineProfile._ymax
protected

Definition at line 133 of file maskStreaks.py.

◆ _ymesh

lsst.meas.algorithms.maskStreaks.LineProfile._ymesh
protected

Definition at line 138 of file maskStreaks.py.

◆ data

lsst.meas.algorithms.maskStreaks.LineProfile.data

Definition at line 131 of file maskStreaks.py.

◆ lineMask

lsst.meas.algorithms.maskStreaks.LineProfile.lineMask

Definition at line 157 of file maskStreaks.py.

◆ lineMaskSize

lsst.meas.algorithms.maskStreaks.LineProfile.lineMaskSize

Definition at line 161 of file maskStreaks.py.

◆ mask

lsst.meas.algorithms.maskStreaks.LineProfile.mask

Definition at line 139 of file maskStreaks.py.

◆ weights

lsst.meas.algorithms.maskStreaks.LineProfile.weights

Definition at line 132 of file maskStreaks.py.


The documentation for this class was generated from the following file: