LSST Applications g0f08755f38+c89d42e150,g1635faa6d4+b6cf076a36,g1653933729+a8ce1bb630,g1a0ca8cf93+4c08b13bf7,g28da252d5a+f33f8200ef,g29321ee8c0+0187be18b1,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+e740673f1a,g5fbc88fb19+17cd334064,g7642f7d749+c89d42e150,g781aacb6e4+a8ce1bb630,g80478fca09+f8b2ab54e1,g82479be7b0+e2bd23ab8b,g858d7b2824+c89d42e150,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+065360aec4,gacf8899fa4+9553554aa7,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gbd46683f8f+ac57cbb13d,gc28159a63d+9634bc57db,gcf0d15dbbd+e37acf7834,gda3e153d99+c89d42e150,gda6a2b7d83+e37acf7834,gdaeeff99f8+1711a396fd,ge2409df99d+cb1e6652d6,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+02b11634a5,w.2024.45
LSST Data Management Base Package
Loading...
Searching...
No Matches
Classes | Functions | Variables
lsst.scarlet.lite.models.parametric Namespace Reference

Classes

class  CartesianFrame
 
class  EllipseFrame
 
class  EllipticalParametricComponent
 
class  ParametricComponent
 

Functions

np.ndarray gaussian2d (np.ndarray params, EllipseFrame ellipse)
 
np.ndarray grad_gaussian2 (np.ndarray input_grad, np.ndarray params, np.ndarray morph, np.ndarray spectrum, EllipseFrame ellipse)
 
np.ndarray circular_gaussian (Sequence[int] center, CartesianFrame frame, float sigma)
 
np.ndarray grad_circular_gaussian (np.ndarray input_grad, np.ndarray params, np.ndarray morph, np.ndarray spectrum, CartesianFrame frame, float sigma)
 
 integrated_gaussian (np.ndarray params, CartesianFrame frame)
 
np.ndarray grad_integrated_gaussian (np.ndarray input_grad, np.ndarray params, np.ndarray morph, np.ndarray spectrum, CartesianFrame frame)
 
np.ndarray bounded_prox (np.ndarray params, np.ndarray proxmin, np.ndarray proxmax)
 
np.ndarray sersic (np.ndarray params, EllipseFrame ellipse)
 
np.ndarray grad_sersic (np.ndarray input_grad, np.ndarray params, np.ndarray morph, np.ndarray spectrum, EllipseFrame ellipse)
 

Variables

int MIN_RADIUS = 1e-20
 
 SQRT_PI_2 = np.sqrt(np.pi / 2)
 
 SERSIC_B1 = gamma.ppf(0.5, 2)
 

Function Documentation

◆ bounded_prox()

np.ndarray lsst.scarlet.lite.models.parametric.bounded_prox ( np.ndarray params,
np.ndarray proxmin,
np.ndarray proxmax )
A bounded proximal operator

This function updates `params` in place.

Parameters
----------
params:
    The array of parameters to constrain.
proxmin:
    The array of minimum values for each parameter.
proxmax:
    The array of maximum values for each parameter.

Returns
-------
result:
    The updated parameters.

Definition at line 535 of file parametric.py.

535def bounded_prox(params: np.ndarray, proxmin: np.ndarray, proxmax: np.ndarray) -> np.ndarray:
536 """A bounded proximal operator
537
538 This function updates `params` in place.
539
540 Parameters
541 ----------
542 params:
543 The array of parameters to constrain.
544 proxmin:
545 The array of minimum values for each parameter.
546 proxmax:
547 The array of maximum values for each parameter.
548
549 Returns
550 -------
551 result:
552 The updated parameters.
553 """
554 cuts = params < proxmin
555 params[cuts] = proxmin[cuts]
556 cuts = params > proxmax
557 params[cuts] = proxmax[cuts]
558 return params
559
560

◆ circular_gaussian()

np.ndarray lsst.scarlet.lite.models.parametric.circular_gaussian ( Sequence[int] center,
CartesianFrame frame,
float sigma )
Model of a circularly symmetric Gaussian

Parameters
----------
center:
    The center of the Gaussian.
frame:
    The frame in which to generate the image of the circular Gaussian
sigma:
    The standard deviation.

Returns
-------
result:
    The image of the circular Gaussian.

Definition at line 388 of file parametric.py.

388def circular_gaussian(center: Sequence[int], frame: CartesianFrame, sigma: float) -> np.ndarray:
389 """Model of a circularly symmetric Gaussian
390
391 Parameters
392 ----------
393 center:
394 The center of the Gaussian.
395 frame:
396 The frame in which to generate the image of the circular Gaussian
397 sigma:
398 The standard deviation.
399
400 Returns
401 -------
402 result:
403 The image of the circular Gaussian.
404 """
405 y0, x0 = center[:2]
406 two_sigma = 2 * sigma
407 r2 = ((frame.x_grid - x0) / two_sigma) ** 2 + ((frame.y_grid - y0) / two_sigma) ** 2
408 return np.exp(-r2)
409
410

◆ gaussian2d()

np.ndarray lsst.scarlet.lite.models.parametric.gaussian2d ( np.ndarray params,
EllipseFrame ellipse )
Model of a 2D elliptical gaussian

Parameters
----------
params:
    The parameters of the function.
    In this case there are none outside of the ellipticity
ellipse:
    The ellipse parameters to scale the radius in all directions.

Returns
-------
result:
    The 2D guassian for the given ellipse parameters

Definition at line 335 of file parametric.py.

335def gaussian2d(params: np.ndarray, ellipse: EllipseFrame) -> np.ndarray:
336 """Model of a 2D elliptical gaussian
337
338 Parameters
339 ----------
340 params:
341 The parameters of the function.
342 In this case there are none outside of the ellipticity
343 ellipse:
344 The ellipse parameters to scale the radius in all directions.
345
346 Returns
347 -------
348 result:
349 The 2D guassian for the given ellipse parameters
350 """
351 return np.exp(-ellipse.r2_grid)
352
353

◆ grad_circular_gaussian()

np.ndarray lsst.scarlet.lite.models.parametric.grad_circular_gaussian ( np.ndarray input_grad,
np.ndarray params,
np.ndarray morph,
np.ndarray spectrum,
CartesianFrame frame,
float sigma )
Gradient of the the component model wrt the Gaussian
morphology parameters

Parameters
----------
input_grad:
    Gradient of the likelihood wrt the component model
params:
    The parameters of the morphology.
morph:
    The model of the morphology.
spectrum:
    The model of the spectrum.
frame:
    The frame in which to generate the image of the circular Gaussian.
sigma:
    The standard deviation.

Definition at line 411 of file parametric.py.

418) -> np.ndarray:
419 """Gradient of the the component model wrt the Gaussian
420 morphology parameters
421
422 Parameters
423 ----------
424 input_grad:
425 Gradient of the likelihood wrt the component model
426 params:
427 The parameters of the morphology.
428 morph:
429 The model of the morphology.
430 spectrum:
431 The model of the spectrum.
432 frame:
433 The frame in which to generate the image of the circular Gaussian.
434 sigma:
435 The standard deviation.
436 """
437 # Calculate the gradient of the likelihod
438 # wrt the Gaussian e^-r**2
439
440 _grad = -morph * np.einsum("i,i...", spectrum, input_grad)
441
442 y0, x0 = params[:2]
443 d_y0 = -2 * np.sum((frame.y_grid - y0) * _grad)
444 d_x0 = -2 * np.sum((frame.x_grid - x0) * _grad)
445 return np.array([d_y0, d_x0], dtype=params.dtype)
446
447

◆ grad_gaussian2()

np.ndarray lsst.scarlet.lite.models.parametric.grad_gaussian2 ( np.ndarray input_grad,
np.ndarray params,
np.ndarray morph,
np.ndarray spectrum,
EllipseFrame ellipse )
Gradient of the the component model wrt the Gaussian
morphology parameters

Parameters
----------
input_grad:
    Gradient of the likelihood wrt the component model
params:
    The parameters of the morphology.
morph:
    The model of the morphology.
spectrum:
    The model of the spectrum.
ellipse:
    The ellipse parameters to scale the radius in all directions.

Definition at line 354 of file parametric.py.

360) -> np.ndarray:
361 """Gradient of the the component model wrt the Gaussian
362 morphology parameters
363
364 Parameters
365 ----------
366 input_grad:
367 Gradient of the likelihood wrt the component model
368 params:
369 The parameters of the morphology.
370 morph:
371 The model of the morphology.
372 spectrum:
373 The model of the spectrum.
374 ellipse:
375 The ellipse parameters to scale the radius in all directions.
376 """
377 # Calculate the gradient of the likelihod
378 # wrt the Gaussian e^-r**2
379 _grad = -morph * np.einsum("i,i...", spectrum, input_grad)
380 d_y0 = ellipse.grad_y0(_grad, True)
381 d_x0 = ellipse.grad_x0(_grad, True)
382 d_sigma_y = ellipse.grad_major(_grad, True)
383 d_sigma_x = ellipse.grad_minor(_grad, True)
384 d_theta = ellipse.grad_theta(_grad, True)
385 return np.array([d_y0, d_x0, d_sigma_y, d_sigma_x, d_theta], dtype=params.dtype)
386
387

◆ grad_integrated_gaussian()

np.ndarray lsst.scarlet.lite.models.parametric.grad_integrated_gaussian ( np.ndarray input_grad,
np.ndarray params,
np.ndarray morph,
np.ndarray spectrum,
CartesianFrame frame )
Gradient of the component model wrt the Gaussian
morphology parameters

Parameters
----------
input_grad:
    Gradient of the likelihood wrt the component model parameters.
params:
    The parameters of the morphology.
morph:
    The model of the morphology.
spectrum:
    The model of the spectrum.
frame:
    The frame in which to generate the image of the circular Gaussian.

Returns
-------
result:
    The gradient of the component morphology.

Definition at line 478 of file parametric.py.

484) -> np.ndarray:
485 """Gradient of the component model wrt the Gaussian
486 morphology parameters
487
488 Parameters
489 ----------
490 input_grad:
491 Gradient of the likelihood wrt the component model parameters.
492 params:
493 The parameters of the morphology.
494 morph:
495 The model of the morphology.
496 spectrum:
497 The model of the spectrum.
498 frame:
499 The frame in which to generate the image of the circular Gaussian.
500
501 Returns
502 -------
503 result:
504 The gradient of the component morphology.
505 """
506 # Calculate the gradient of the likelihood
507 # wrt the Gaussian e^-r**2
508 _grad = np.einsum("i,i...", spectrum, input_grad)
509
510 # Extract the parameters
511 y0, x0, sigma = params
512 # define useful constants
513 x = frame.x_grid - x0
514 y = frame.y_grid - y0
515 c = 0.5 / sigma**2
516 sqrt_c = np.sqrt(c)
517 # Add a small constant to the radius to prevent a divergence at r==0
518 r = np.sqrt(x**2 + y**2) + MIN_RADIUS
519 # Shift half a pixel in each direction for the integration
520 r1 = r - 0.5
521 r2 = r + 0.5
522 # Calculate the gradient of the ERF wrt. each shifted radius
523 d_model1 = np.exp(-c * r1**2)
524 d_model2 = np.exp(-c * r2**2)
525 # Calculate the gradients of the parameters
526 d_x0 = np.sum(-x / r * (d_model2 - d_model1) * _grad)
527 d_y0 = np.sum(-y / r * (d_model2 - d_model1) * _grad)
528 d_sigma1 = -(r1 * d_model1 / sigma - SQRT_PI_2 * erf(r1 * sqrt_c))
529 d_sigma2 = -(r2 * d_model2 / sigma - SQRT_PI_2 * erf(r2 * sqrt_c))
530 d_sigma = np.sum((d_sigma2 - d_sigma1) * _grad)
531
532 return np.array([d_y0, d_x0, d_sigma])
533
534

◆ grad_sersic()

np.ndarray lsst.scarlet.lite.models.parametric.grad_sersic ( np.ndarray input_grad,
np.ndarray params,
np.ndarray morph,
np.ndarray spectrum,
EllipseFrame ellipse )
Gradient of the component model wrt the morphology parameters

Parameters
----------
input_grad:
    Gradient of the likelihood wrt the component model
params:
    The parameters of the morphology.
morph:
    The model of the morphology.
spectrum:
    The model of the spectrum.
ellipse:
    The ellipse parameters to scale the radius in all directions.

Definition at line 589 of file parametric.py.

595) -> np.ndarray:
596 """Gradient of the component model wrt the morphology parameters
597
598 Parameters
599 ----------
600 input_grad:
601 Gradient of the likelihood wrt the component model
602 params:
603 The parameters of the morphology.
604 morph:
605 The model of the morphology.
606 spectrum:
607 The model of the spectrum.
608 ellipse:
609 The ellipse parameters to scale the radius in all directions.
610 """
611 n = params[5]
612 bn = gamma.ppf(0.5, 2 * n)
613 if n == 1:
614 # Use a simplified model for faster calculation
615 d_exp = -SERSIC_B1 * morph
616 else:
617 r = ellipse.r_grid
618 d_exp = -bn / n * morph * r ** (1 / n - 1)
619
620 _grad = np.einsum("i,i...", spectrum, input_grad)
621 d_n = np.sum(_grad * bn * morph * ellipse.r_grid ** (1 / n) * np.log10(ellipse.r_grid) / n**2)
622 _grad = _grad * d_exp
623 d_y0 = ellipse.grad_y0(_grad, False)
624 d_x0 = ellipse.grad_x0(_grad, False)
625 d_sigma_y = ellipse.grad_major(_grad, False)
626 d_sigma_x = ellipse.grad_minor(_grad, False)
627 d_theta = ellipse.grad_theta(_grad, False)
628 return np.array([d_y0, d_x0, d_sigma_y, d_sigma_x, d_theta, d_n], dtype=params.dtype)
629
630

◆ integrated_gaussian()

lsst.scarlet.lite.models.parametric.integrated_gaussian ( np.ndarray params,
CartesianFrame frame )
Model of a circularly symmetric Gaussian integrated over pixels

This differs from `circularGaussian` because the gaussian function
is integrated over each pixel to replicate the pixelated image
version of a Gaussian function.

Parameters
----------
params:
    The center of the Gaussian.
frame:
    The frame in which to generate the image of the circular Gaussian

Returns
-------
result:
    The image of the circular Gaussian.

Definition at line 448 of file parametric.py.

448def integrated_gaussian(params: np.ndarray, frame: CartesianFrame):
449 """Model of a circularly symmetric Gaussian integrated over pixels
450
451 This differs from `circularGaussian` because the gaussian function
452 is integrated over each pixel to replicate the pixelated image
453 version of a Gaussian function.
454
455 Parameters
456 ----------
457 params:
458 The center of the Gaussian.
459 frame:
460 The frame in which to generate the image of the circular Gaussian
461
462 Returns
463 -------
464 result:
465 The image of the circular Gaussian.
466 """
467 # Unpack the parameters and define constants
468 y0, x0, sigma = params
469 r = np.sqrt((frame.x_grid - x0) ** 2 + (frame.y_grid - y0) ** 2)
470 sqrt_c = 1 / np.sqrt(2) / sigma
471 # Integrate from half a pixel left and right
472 lhs = erf((r - 0.5) * sqrt_c)
473 rhs = erf((r + 0.5) * sqrt_c)
474 z = 0.5 * np.sqrt(np.pi) / sqrt_c * (rhs - lhs)
475 return z
476
477

◆ sersic()

np.ndarray lsst.scarlet.lite.models.parametric.sersic ( np.ndarray params,
EllipseFrame ellipse )
Generate a Sersic Model.

Parameters
----------
params:
    The parameters of the function.
    In this case the only parameter is the sersic index ``n``.
ellipse:
    The ellipse parameters to scale the radius in all directions.

Returns
-------
result:
    The model for the given ellipse parameters

Definition at line 561 of file parametric.py.

561def sersic(params: np.ndarray, ellipse: EllipseFrame) -> np.ndarray:
562 """Generate a Sersic Model.
563
564 Parameters
565 ----------
566 params:
567 The parameters of the function.
568 In this case the only parameter is the sersic index ``n``.
569 ellipse:
570 The ellipse parameters to scale the radius in all directions.
571
572 Returns
573 -------
574 result:
575 The model for the given ellipse parameters
576 """
577 (n,) = params
578
579 r = ellipse.r_grid
580
581 if n == 1:
582 result = np.exp(-SERSIC_B1 * r)
583 else:
584 bn = gamma.ppf(0.5, 2 * n)
585 result = np.exp(-bn * (r ** (1 / n) - 1))
586 return result
587
588

Variable Documentation

◆ MIN_RADIUS

int lsst.scarlet.lite.models.parametric.MIN_RADIUS = 1e-20

Definition at line 51 of file parametric.py.

◆ SERSIC_B1

lsst.scarlet.lite.models.parametric.SERSIC_B1 = gamma.ppf(0.5, 2)

Definition at line 57 of file parametric.py.

◆ SQRT_PI_2

lsst.scarlet.lite.models.parametric.SQRT_PI_2 = np.sqrt(np.pi / 2)

Definition at line 54 of file parametric.py.