LSST Applications g00d0e8bbd7+8c5ae1fdc5,g013ef56533+603670b062,g083dd6704c+2e189452a7,g199a45376c+0ba108daf9,g1c5cce2383+bc9f6103a4,g1fd858c14a+cd69ed4fc1,g210f2d0738+c4742f2e9e,g262e1987ae+612fa42d85,g29ae962dfc+83d129e820,g2cef7863aa+aef1011c0b,g35bb328faa+8c5ae1fdc5,g3fd5ace14f+5eaa884f2a,g47891489e3+e32160a944,g53246c7159+8c5ae1fdc5,g5b326b94bb+dcc56af22d,g64539dfbff+c4742f2e9e,g67b6fd64d1+e32160a944,g74acd417e5+c122e1277d,g786e29fd12+668abc6043,g87389fa792+8856018cbb,g88cb488625+47d24e4084,g89139ef638+e32160a944,g8d7436a09f+d14b4ff40a,g8ea07a8fe4+b212507b11,g90f42f885a+e1755607f3,g97be763408+34be90ab8c,g98df359435+ec1fa61bf1,ga2180abaac+8c5ae1fdc5,ga9e74d7ce9+43ac651df0,gbf99507273+8c5ae1fdc5,gc2a301910b+c4742f2e9e,gca7fc764a6+e32160a944,gd7ef33dd92+e32160a944,gdab6d2f7ff+c122e1277d,gdb1e2cdc75+1b18322db8,ge410e46f29+e32160a944,ge41e95a9f2+c4742f2e9e,geaed405ab2+0d91c11c6d,w.2025.44
LSST Data Management Base Package
Loading...
Searching...
No Matches
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 540 of file parametric.py.

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

◆ 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 393 of file parametric.py.

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

◆ 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 340 of file parametric.py.

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

◆ 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 416 of file parametric.py.

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

◆ 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 359 of file parametric.py.

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

◆ 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 483 of file parametric.py.

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

◆ 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 594 of file parametric.py.

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

◆ 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 453 of file parametric.py.

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

◆ 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 566 of file parametric.py.

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

Variable Documentation

◆ MIN_RADIUS

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

Definition at line 56 of file parametric.py.

◆ SERSIC_B1

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

Definition at line 62 of file parametric.py.

◆ SQRT_PI_2

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

Definition at line 59 of file parametric.py.