LSST Applications g0265f82a02+d6b5cd48b5,g02d81e74bb+80768bd682,g04242d3e92+8eaa23c173,g06b2ea86fd+734f9505a2,g2079a07aa2+14824f138e,g212a7c68fe+5f4fc2ea00,g2305ad1205+293ab1327e,g2bbee38e9b+d6b5cd48b5,g337abbeb29+d6b5cd48b5,g3ddfee87b4+8eaa23c173,g487adcacf7+abec5a19c5,g50ff169b8f+5929b3527e,g52b1c1532d+a6fc98d2e7,g591dd9f2cf+97ef3b4495,g5a732f18d5+66d966b544,g5d7b63bc56+636c3c3fd8,g64a986408d+80768bd682,g858d7b2824+80768bd682,g8a8a8dda67+a6fc98d2e7,g99cad8db69+6282a5f541,g9ddcbc5298+d4bad12328,ga1e77700b3+246acaaf9c,ga8c6da7877+9e3c062e8e,gb0e22166c9+3863383f4c,gb6a65358fc+d6b5cd48b5,gba4ed39666+9664299f35,gbb8dafda3b+60f904e7bc,gc120e1dc64+1bf26d0180,gc28159a63d+d6b5cd48b5,gcf0d15dbbd+8eaa23c173,gd2a12a3803+f8351bc914,gdaeeff99f8+a38ce5ea23,ge79ae78c31+d6b5cd48b5,gee10cc3b42+a6fc98d2e7,gf1cff7945b+80768bd682,v24.1.5.rc1
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
lsst.pipe.tasks.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 189 of file maskStreaks.py.

Constructor & Destructor Documentation

◆ __init__()

lsst.pipe.tasks.maskStreaks.LineProfile.__init__ ( self,
data,
weights,
line = None )

Definition at line 213 of file maskStreaks.py.

213 def __init__(self, data, weights, line=None):
214 self.data = data
215 self.weights = weights
216 self._ymax, self._xmax = data.shape
217 self._dtype = data.dtype
218 xrange = np.arange(self._xmax) - self._xmax / 2.
219 yrange = np.arange(self._ymax) - self._ymax / 2.
220 self._rhoMax = ((0.5 * self._ymax)**2 + (0.5 * self._xmax)**2)**0.5
221 self._xmesh, self._ymesh = np.meshgrid(xrange, yrange)
222 self.mask = (weights != 0)
223
224 self._initLine = line
225 self.setLineMask(line)
226

Member Function Documentation

◆ _lineChi2()

lsst.pipe.tasks.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 328 of file maskStreaks.py.

328 def _lineChi2(self, line, grad=True):
329 """Construct the chi2 between the data and the model.
330
331 Parameters
332 ----------
333 line : `Line`
334 `Line` parameters for which to build model and calculate chi2.
335 grad : `bool`, optional
336 Whether or not to return the gradient and hessian.
337
338 Returns
339 -------
340 reducedChi : `float`
341 Reduced chi2 of the model.
342 reducedDChi : `np.ndarray`
343 Derivative of the chi2 with respect to rho, theta, invSigma.
344 reducedHessianChi : `np.ndarray`
345 Hessian of the chi2 with respect to rho, theta, invSigma.
346 """
347 # Calculate chi2
348 model, dModel = self._makeMaskedProfile(line)
349 chi2 = (self._maskWeights * (self._maskData - model)**2).sum()
350 if not grad:
351 return chi2.sum() / self.lineMaskSize
352
353 # Calculate derivative and Hessian of chi2
354 derivChi2 = ((-2 * self._maskWeights * (self._maskData - model))[None, :] * dModel).sum(axis=1)
355 hessianChi2 = (2 * self._maskWeights * dModel[:, None, :] * dModel[None, :, :]).sum(axis=2)
356
357 reducedChi = chi2 / self.lineMaskSize
358 reducedDChi = derivChi2 / self.lineMaskSize
359 reducedHessianChi = hessianChi2 / self.lineMaskSize
360 return reducedChi, reducedDChi, reducedHessianChi
361

◆ _makeMaskedProfile()

lsst.pipe.tasks.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 250 of file maskStreaks.py.

250 def _makeMaskedProfile(self, line, fitFlux=True):
251 """Construct the line model in the masked region and calculate its
252 derivatives.
253
254 Parameters
255 ----------
256 line : `Line`
257 Parameters of line profile for which to make profile in the masked
258 region.
259 fitFlux : `bool`
260 Fit the amplitude of the line profile to the data.
261
262 Returns
263 -------
264 model : `np.ndarray`
265 Model in the masked region.
266 dModel : `np.ndarray`
267 Derivative of the model in the masked region.
268 """
269 invSigma = line.sigma**-1
270 # Calculate distance between pixels and line
271 radtheta = np.deg2rad(line.theta)
272 costheta = np.cos(radtheta)
273 sintheta = np.sin(radtheta)
274 distance = (costheta * self._mxmesh + sintheta * self._mymesh - line.rho)
275 distanceSquared = distance**2
276
277 # Calculate partial derivatives of distance
278 drad = np.pi / 180
279 dDistanceSqdRho = 2 * distance * (-np.ones_like(self._mxmesh))
280 dDistanceSqdTheta = (2 * distance * (-sintheta * self._mxmesh + costheta * self._mymesh) * drad)
281
282 # Use pixel-line distances to make Moffat profile
283 profile = (1 + distanceSquared * invSigma**2)**-2.5
284 dProfile = -2.5 * (1 + distanceSquared * invSigma**2)**-3.5
285
286 if fitFlux:
287 # Calculate line flux from profile and data
288 flux = ((self._maskWeights * self._maskData * profile).sum()
289 / (self._maskWeights * profile**2).sum())
290 else:
291 # Approximately normalize the line
292 flux = invSigma**-1
293 if np.isnan(flux):
294 flux = 0
295
296 model = flux * profile
297
298 # Calculate model derivatives
299 fluxdProfile = flux * dProfile
300 fluxdProfileInvSigma = fluxdProfile * invSigma**2
301 dModeldRho = fluxdProfileInvSigma * dDistanceSqdRho
302 dModeldTheta = fluxdProfileInvSigma * dDistanceSqdTheta
303 dModeldInvSigma = fluxdProfile * distanceSquared * 2 * invSigma
304
305 dModel = np.array([dModeldRho, dModeldTheta, dModeldInvSigma])
306 return model, dModel
307

◆ fit()

lsst.pipe.tasks.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 362 of file maskStreaks.py.

362 def fit(self, dChi2Tol=0.1, maxIter=100, log=None):
363 """Perform Newton-Raphson minimization to find line parameters.
364
365 This method takes advantage of having known derivative and Hessian of
366 the multivariate function to quickly and efficiently find the minimum.
367 This is more efficient than the scipy implementation of the Newton-
368 Raphson method, which doesn't take advantage of the Hessian matrix. The
369 method here also performs a line search in the direction of the steepest
370 derivative at each iteration, which reduces the number of iterations
371 needed.
372
373 Parameters
374 ----------
375 dChi2Tol : `float`, optional
376 Change in Chi2 tolerated for fit convergence.
377 maxIter : `int`, optional
378 Maximum number of fit iterations allowed. The fit should converge in
379 ~10 iterations, depending on the value of dChi2Tol, but this
380 maximum provides a backup.
381 log : `lsst.utils.logging.LsstLogAdapter`, optional
382 Logger to use for reporting more details for fitting failures.
383
384 Returns
385 -------
386 outline : `np.ndarray`
387 Coordinates and inverse width of fit line.
388 chi2 : `float`
389 Reduced Chi2 of model fit to data.
390 fitFailure : `bool`
391 Boolean where `False` corresponds to a successful fit.
392 """
393 # Do minimization on inverse of sigma to simplify derivatives:
394 x = np.array([self._initLine.rho, self._initLine.theta, self._initLine.sigma**-1])
395
396 dChi2 = 1
397 iter = 0
398 oldChi2 = 0
399 fitFailure = False
400
401 def line_search(c, dx):
402 testx = x - c * dx
403 testLine = Line(testx[0], testx[1], testx[2]**-1)
404 return self._lineChi2(testLine, grad=False)
405
406 while abs(dChi2) > dChi2Tol:
407 line = Line(x[0], x[1], x[2]**-1)
408 chi2, b, A = self._lineChi2(line)
409 if chi2 == 0:
410 break
411 if not np.isfinite(A).all():
412 fitFailure = True
413 if log is not None:
414 log.warning("Hessian matrix has non-finite elements.")
415 break
416 dChi2 = oldChi2 - chi2
417 try:
418 cholesky = scipy.linalg.cho_factor(A)
419 except np.linalg.LinAlgError:
420 fitFailure = True
421 if log is not None:
422 log.warning("Hessian matrix is not invertible.")
423 break
424 dx = scipy.linalg.cho_solve(cholesky, b)
425
426 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05)
427 x -= factor * dx
428 if (abs(x[0]) > 1.5 * self._rhoMax) or (iter > maxIter):
429 fitFailure = True
430 break
431 oldChi2 = chi2
432 iter += 1
433
434 outline = Line(x[0], x[1], abs(x[2])**-1)
435
436 return outline, chi2, fitFailure
437
438

◆ makeProfile()

lsst.pipe.tasks.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 308 of file maskStreaks.py.

308 def makeProfile(self, line, fitFlux=True):
309 """Construct the line profile model.
310
311 Parameters
312 ----------
313 line : `Line`
314 Parameters of the line profile to model.
315 fitFlux : `bool`, optional
316 Fit the amplitude of the line profile to the data.
317
318 Returns
319 -------
320 finalModel : `np.ndarray`
321 Model for line profile.
322 """
323 model, _ = self._makeMaskedProfile(line, fitFlux=fitFlux)
324 finalModel = np.zeros((self._ymax, self._xmax), dtype=self._dtype)
325 finalModel[self.lineMask] = model
326 return finalModel
327

◆ setLineMask()

lsst.pipe.tasks.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 227 of file maskStreaks.py.

227 def setLineMask(self, line):
228 """Set mask around the image region near the line.
229
230 Parameters
231 ----------
232 line : `Line`
233 Parameters of line in the image.
234 """
235 if line:
236 # Only fit pixels within 5 sigma of the estimated line
237 radtheta = np.deg2rad(line.theta)
238 distance = (np.cos(radtheta) * self._xmesh + np.sin(radtheta) * self._ymesh - line.rho)
239 m = (abs(distance) < 5 * line.sigma)
240 self.lineMask = self.mask & m
241 else:
242 self.lineMask = np.copy(self.mask)
243
244 self.lineMaskSize = self.lineMask.sum()
245 self._maskData = self.data[self.lineMask]
246 self._maskWeights = self.weights[self.lineMask]
247 self._mxmesh = self._xmesh[self.lineMask]
248 self._mymesh = self._ymesh[self.lineMask]
249

Member Data Documentation

◆ _dtype

lsst.pipe.tasks.maskStreaks.LineProfile._dtype
protected

Definition at line 217 of file maskStreaks.py.

◆ _initLine

lsst.pipe.tasks.maskStreaks.LineProfile._initLine
protected

Definition at line 224 of file maskStreaks.py.

◆ _maskData

lsst.pipe.tasks.maskStreaks.LineProfile._maskData
protected

Definition at line 245 of file maskStreaks.py.

◆ _maskWeights

lsst.pipe.tasks.maskStreaks.LineProfile._maskWeights
protected

Definition at line 246 of file maskStreaks.py.

◆ _mxmesh

lsst.pipe.tasks.maskStreaks.LineProfile._mxmesh
protected

Definition at line 247 of file maskStreaks.py.

◆ _mymesh

lsst.pipe.tasks.maskStreaks.LineProfile._mymesh
protected

Definition at line 248 of file maskStreaks.py.

◆ _rhoMax

lsst.pipe.tasks.maskStreaks.LineProfile._rhoMax
protected

Definition at line 220 of file maskStreaks.py.

◆ _xmax

lsst.pipe.tasks.maskStreaks.LineProfile._xmax
protected

Definition at line 216 of file maskStreaks.py.

◆ _xmesh

lsst.pipe.tasks.maskStreaks.LineProfile._xmesh
protected

Definition at line 221 of file maskStreaks.py.

◆ _ymax

lsst.pipe.tasks.maskStreaks.LineProfile._ymax
protected

Definition at line 216 of file maskStreaks.py.

◆ _ymesh

lsst.pipe.tasks.maskStreaks.LineProfile._ymesh
protected

Definition at line 221 of file maskStreaks.py.

◆ data

lsst.pipe.tasks.maskStreaks.LineProfile.data

Definition at line 214 of file maskStreaks.py.

◆ lineMask

lsst.pipe.tasks.maskStreaks.LineProfile.lineMask

Definition at line 240 of file maskStreaks.py.

◆ lineMaskSize

lsst.pipe.tasks.maskStreaks.LineProfile.lineMaskSize

Definition at line 244 of file maskStreaks.py.

◆ mask

lsst.pipe.tasks.maskStreaks.LineProfile.mask

Definition at line 222 of file maskStreaks.py.

◆ weights

lsst.pipe.tasks.maskStreaks.LineProfile.weights

Definition at line 215 of file maskStreaks.py.


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