socca.models¶
- class socca.models.Component(**kwargs)[source]¶
Bases:
objectBase class for all model components.
This abstract base class provides the infrastructure for model components, including parameter management, unique identification, and formatted output. All specific profile types (Beta, Sersic, Point, etc.) inherit from this class.
- Variables:
idcls (int (class variable)) – Class-level counter for assigning unique IDs to component instances.
id (str) – Unique identifier for this component instance (e.g., ‘comp_00’, ‘comp_01’).
positive (bool) – Whether to enforce positivity constraint on this component.
hyper (list of str) – Names of hyperparameters (not fitted, control behavior).
units (dict) – Mapping from parameter names to their physical units.
description (dict) – Mapping from parameter names to their descriptions.
- idcls = 0¶
- __init__(**kwargs)[source]¶
Initialize a component with unique ID and parameter tracking structures.
- Parameters:
**kwargs (dict) –
Keyword arguments for component initialization.
- positivebool, optional
Override default positivity constraint.
Notes
The component ID is automatically assigned and incremented for each new instance. This ensures unique identification when multiple components are combined in a model.
- parameters()[source]¶
Print a formatted table of component parameters and hyperparameters.
Displays all parameters (excluding reserved keys) with their values, units, and descriptions in a human-readable format. Separates regular parameters from hyperparameters.
Notes
Output format:
Model parameters¶
parameter_name [units] : value | description
Hyperparameters¶
hyperparam_name [units] : value | description
Parameters are model quantities to be fitted or fixed
Hyperparameters control computation but are not fitted
Values are shown in scientific notation (e.g., 1.5000E-02)
None values are displayed as “None”
Examples
>>> from socca.models import Beta >>> beta = Beta(xc=180.5, yc=45.2, rc=0.01, Ic=100, beta=0.5) >>> beta.parameters() Model parameters ================ xc [deg] : 1.8050E+02 | Right ascension of centroid yc [deg] : 4.5200E+01 | Declination of centroid rc [deg] : 1.0000E-02 | Core radius Ic [image] : 1.0000E+02 | Central surface brightness beta [] : 5.0000E-01 | Slope parameter theta [rad] : 0.0000E+00 | Position angle (east from north) e [] : 0.0000E+00 | Projected axis ratio cbox [] : 0.0000E+00 | Projected boxiness
- parlist()[source]¶
Return a list of parameter names for this component.
- Returns:
Names of all model parameters, as defined in the units dict.
- Return type:
list of str
Notes
This method is used internally by Model.addcomponent() to register parameters when adding a component to a composite model.
Examples
>>> from socca.models import Beta >>> beta = Beta(xc=180.5, yc=45.2) >>> beta.parlist() ['xc', 'yc', 'theta', 'e', 'cbox', 'rc', 'Ic', 'beta']
- addparameter(name, value=None, units='', description='')[source]¶
Add a new parameter to the component with metadata.
- Parameters:
name (str) – Name of the parameter to add.
value (float, int, Distribution, or callable, optional) – Parameter value. Can be a fixed value, a numpyro prior distribution, or a function for tied parameters. Default is None.
units (str, optional) – Physical units of the parameter (e.g., ‘deg’, ‘rad’, ‘image’). Default is empty string.
description (str, optional) – Human-readable description of the parameter. Default is empty string.
Notes
This method dynamically adds a parameter attribute to the component instance and registers its metadata (units and description).
Examples
>>> from socca.models import Component >>> comp = Component() >>> comp.addparameter('amplitude', value=1.0, units='Jy', description='Source flux') >>> comp.amplitude 1.0 >>> comp.units['amplitude'] 'Jy'
- removeparameter(name)[source]¶
Remove a parameter from the component.
- Parameters:
name (str) – Name of the parameter to remove.
- Raises:
AttributeError – If the parameter does not exist on this component.
Examples
>>> from socca.models import Beta >>> beta = Beta(xc=180.5, yc=45.2) >>> beta.addparameter('custom', value=1.0, units='Jy') >>> beta.removeparameter('custom')
- class socca.models.Model(prof=None, positive=None)[source]¶
Bases:
objectComposite model container for combining multiple profile components.
The Model class provides a flexible framework for constructing complex astronomical models by combining individual profile components (Beta, Sersic, Point, etc.). It manages parameter priors, positivity constraints, and model evaluation.
- __init__(prof=None, positive=None)[source]¶
Initialize a composite model with optional initial profile component.
- Parameters:
prof (Profile or Component, optional) – Initial profile component to add to the model. Can be any profile instance (Beta, Sersic, Point, etc.). If None, creates an empty model to which components can be added later.
positive (bool, optional) – Whether to enforce positivity constraint on the profile. If None, uses the default positivity setting from the profile itself.
- Variables:
ncomp (int) – Number of components in the model.
priors (dict) – Dictionary mapping parameter names to their values or prior distributions.
params (list of str) – List of all parameter names in the model.
paridx (list of int) – Indices of parameters with prior distributions (free parameters).
positive (list of bool) – Positivity constraints for each component.
tied (list of bool) – Indicates which parameters are tied to other parameters.
type (list of str) – Class names of each component in the model.
units (dict) – Dictionary mapping parameter names to their physical units.
Examples
>>> from socca.models import Model, Beta >>> # Create empty model and add component later >>> model = Model() >>> model.addcomponent(Beta()) >>> # Or initialize with a profile directly >>> model = Model(Beta())
- addcomponent(prof, positive=None)[source]¶
Add a profile component to the composite model.
This method incorporates a new profile component into the model, handling parameter registration, prior assignments, and constraint management. Components are indexed sequentially and their parameters are prefixed with the component index.
- Parameters:
- Raises:
ValueError – If any parameter in the profile is set to None without a valid value or prior distribution.
Notes
Parameters are stored with names like ‘comp_00_xc’, ‘comp_01_rc’, etc.
Parameters can be fixed values, prior distributions, or tied (functions).
Only parameters with prior distributions are considered free parameters during fitting.
Tied parameters are evaluated as functions of other parameters.
Examples
>>> from socca.models import Model, Beta, Point >>> model = Model() >>> # Add a Beta profile >>> model.addcomponent(Beta(xc=180.5, yc=45.2, rc=0.01)) >>> # Add a point source >>> model.addcomponent(Point(xc=180.6, yc=45.3), positive=True) >>> print(model.ncomp) 2
- parameters(freeonly=False)[source]¶
Print a formatted table of model parameters with their values and units.
This method displays all model parameters, showing fixed values, prior distributions, or tied parameter relationships in a human-readable format.
- Parameters:
freeonly (bool, optional) – If True, only display parameters with prior distributions (free parameters to be fitted). If False (default), display all parameters including fixed values and tied parameters.
Notes
Parameter values are displayed as:
Fixed values: shown in scientific notation (e.g., 1.5000E-02)
Prior distributions: shown as “Distribution: DistributionName”
Tied parameters: shown as “Tied parameter”
The output format is: parameter_name [units] : value
Examples
>>> from socca.models import Model, Beta >>> import numpyro.distributions as dist >>> model = Model(Beta(xc=dist.Uniform(180, 181), yc=45.2, rc=0.01)) >>> model.parameters() Model parameters ================ comp_00_xc [deg] : Distribution: Uniform comp_00_yc [deg] : 4.5200E+01 comp_00_rc [deg] : 1.0000E-02 ... >>> model.parameters(freeonly=True) Model parameters ================ comp_00_xc [deg] : Distribution: Uniform
- getmap(img, convolve=None, addbackground=False)[source]¶
Generate a model image map using fixed parameter values only.
This method evaluates the model on the provided image grid using only fixed parameter values. All parameters must be fixed; parameters with prior distributions will raise an error.
- Parameters:
img (Image) – Image object containing the grid, PSF, and WCS information.
convolve (bool, optional) – If True, return the PSF-convolved model. If False (default), return the unconvolved model. If None, defaults to False.
addbackground (bool, optional) – If True and convolve is True, include the background component in the output. If False (default), exclude background. Note that background is only added to convolved maps.
- Returns:
Model image map on the same grid as img.data. Shape matches img.data.shape.
- Return type:
ndarray
- Raises:
ValueError – If any parameter has a prior distribution instead of a fixed value.
- Warns:
UserWarning – If addbackground=True but convolve=False, warns that background is only added to convolved maps.
Notes
All model parameters must be fixed values (float or int).
Tied parameters are automatically evaluated from their dependencies.
For use during fitting with sampled parameters, use getmodel() instead.
Examples
>>> from socca.models import Model, Beta >>> from socca.data import Image >>> model = Model(Beta(xc=180.5, yc=45.2, rc=0.01, Ic=100, beta=0.5)) >>> img = Image('observation.fits') >>> # Get unconvolved model >>> model_map = model.getmap(img) >>> # Get PSF-convolved model >>> convolved_map = model.getmap(img, convolve=True)
- getmodel(img, pp, doresp=False, doexp=False, component=None)[source]¶
Compute the full model image with PSF convolution, response, and exposure.
This is the main forward modeling function used during fitting. It evaluates all model components, applies PSF convolution, instrumental response, and exposure corrections as specified.
- Parameters:
img (Image) – Image object containing data, grid, PSF, response map, and exposure map.
pp (array_like or dict) – Parameter values. If array_like, contains values for parameters with prior distributions, in the order they appear in paridx. If dict, contains all parameter values with keys matching self.params.
doresp (bool, optional) – If True, multiply by the instrumental response map (img.response). Default is False.
doexp (bool, optional) – If True, multiply by the exposure map (img.exposure). Default is False.
component (None, str, int, list, or Profile, optional) –
Component(s) to include in the model computation. Can be:
None: Include all components (default)
str: Single component name (e.g., ‘comp_00’)
int: Component index (e.g., 0 for first component)
list: Multiple components as names, indices, or Profile objects
Profile: Object with id attribute specifying the component
Default is None (all components).
- Returns:
mraw (ndarray) – Unconvolved model image (raw surface brightness distribution).
msmo (ndarray) – Final model image after PSF convolution, response, and exposure corrections, including background.
mbkg (ndarray) – Background component only, with exposure applied if doexp=True.
mneg (ndarray) – Mask indicating pixels where model would be negative (if positivity constraints are violated). 1.0 where negative, 0.0 otherwise.
Notes
Fixed parameters are taken from self.priors
Free parameters are extracted from pp array
Tied parameters are computed from their functional dependencies
Point sources are handled in Fourier space for efficiency
Background components are not PSF-convolved
Disk components use 3D line-of-sight integration
The positivity mask (mneg) can be used to penalize negative values
Examples
>>> from socca.models import Model, Beta >>> import numpyro.distributions as dist >>> model = Model(Beta(xc=dist.Uniform(180, 181), yc=45.2)) >>> img = Image('observation.fits') >>> # During fitting, pp contains sampled parameter values >>> pp = [180.5] # value for xc >>> mraw, msmo, mbkg, mneg = model.getmodel(img, pp, doresp=True, doexp=True)
- socca.models.zoo()[source]¶
Print available model profiles in the socca library.
This function displays a list of all available profile models that can be used for fitting astronomical images, including various analytical profiles for galaxies, point sources, and background components.
Notes
Available models include:
Beta: Beta profile (power-law density)
gNFW: Generalized Navarro-Frenk-White profile
Power: Power-law profile
TopHat: Uniform top-hat profile
Sersic: Sersic profile for elliptical galaxies
Gaussian: Gaussian profile
Exponential: Exponential disk profile
PolyExponential: Polynomial-exponential profile
PolyExpoRefact: Refactored polynomial-exponential profile
ModExponential: Modified exponential profile
Point: Point source model
Background: Polynomial background model
Disk: 3D disk model with finite thickness
Examples
>>> import socca.models as models >>> models.zoo() Beta gNFW Power TopHat Sersic Gaussian Exponential PolyExponential PolyExpoRefact ModExponential Point Background Disk
- class socca.models.Profile(**kwargs)[source]¶
Bases:
ComponentBase class for 2D surface brightness profiles.
This class extends Component to provide common coordinate transformation and projection parameters (position, orientation, ellipticity, boxiness) for all 2D surface brightness profiles.
- Parameters:
**kwargs (dict) –
Keyword arguments for profile initialization including:
- xcfloat, optional
Right ascension of profile centroid (deg).
- ycfloat, optional
Declination of profile centroid (deg).
- thetafloat, optional
Position angle measured east from north (rad).
- efloat, optional
Ellipticity, defined as 1 - b/a where b/a is axis ratio (0 <= e < 1).
- cboxfloat, optional
Boxiness parameter (0 = elliptical, >0 = boxy, <0 = disky).
- positivebool, optional
Whether to enforce positivity constraint.
- Variables:
yc (xc,) – Centroid coordinates in celestial degrees.
theta (float) – Position angle in radians.
e (float) – Ellipticity parameter.
cbox (float) – Boxiness/diskyness parameter.
Notes
All Profile subclasses must implement the abstract profile() method that defines the radial surface brightness distribution.
- __init__(**kwargs)[source]¶
Initialize a profile with standard coordinate and shape parameters.
- Parameters:
**kwargs (dict) – Keyword arguments including xc, yc, theta, e, cbox, positive. See class docstring for parameter descriptions.
- property truncation¶
- property e¶
- abstractmethod profile(r)[source]¶
Abstract method defining the radial profile function.
Subclasses must implement this method to define the surface brightness distribution as a function of radius and other profile-specific parameters.
- Parameters:
r (ndarray) – Radial coordinate in degrees.
- Returns:
Surface brightness values at each radius.
- Return type:
ndarray
- getmap(img, convolve=False)[source]¶
Generate an image map from the profile on the given grid.
Evaluates the profile on the image grid and optionally convolves with the PSF. All parameters must have fixed values.
- Parameters:
img (Image) – Image object containing grid, PSF, and WCS information.
convolve (bool, optional) – If True, convolve the model with the PSF. Default is False.
- Returns:
2D array of surface brightness values on the image grid. Shape matches img.data.shape.
- Return type:
ndarray
- Raises:
ValueError – If any parameter is a prior distribution or set to None.
- Warns:
UserWarning – If convolve=True but no PSF is defined in img.
Notes
All profile parameters must be fixed values (not distributions).
The elliptical grid is computed using position, ellipticity, and PA.
Convolution is performed in Fourier space for efficiency.
Examples
>>> from socca.models import Beta >>> from socca.data import Image >>> beta = Beta(xc=180.5, yc=45.2, rc=0.01, Ic=100, beta=0.5) >>> img = Image('observation.fits') >>> model_map = beta.getmap(img, convolve=True)
- static getgrid(grid, xc, yc, theta=0.0, e=0.0, cbox=0.0)[source]¶
Compute elliptical radius grid with rotation and projection.
This static JIT-compiled method transforms celestial coordinates to elliptical radius values accounting for position angle, ellipticity, and boxiness. Used internally by profile evaluation.
- Parameters:
grid (Grid) – Grid object with .x and .y celestial coordinate arrays (deg).
xc (float) – Right ascension of centroid (deg).
yc (float) – Declination of centroid (deg).
theta (float, optional) – Position angle east from north (rad). Default is 0.
e (float, optional) – Ellipticity (1 - b/a). Default is 0 (circular).
cbox (float, optional) – Boxiness parameter. Default is 0 (elliptical).
- Returns:
Elliptical radius grid in degrees. Same shape as grid.x and grid.y.
- Return type:
ndarray
Notes
The transformation accounts for:
Spherical geometry (cos(dec) correction)
Position angle rotation
Ellipticity via axis ratio correction
Generalized elliptical isophotes with boxiness
The elliptical radius is computed as: r = [(abs(x’)^(2+c) + abs(y’/(1-e))^(2+c)]^(1/(2+c)) where c is the boxiness parameter.
This function is JIT-compiled for performance during model evaluation.
- refactor()[source]¶
Return a refactored version of the profile with equivalent parameterization.
For most profiles, this returns a copy of the profile with the same parameters. Some profiles (e.g., PolyExpoRefact) override this to convert to an alternative parameterization.
- Returns:
A new profile instance with refactored parameterization.
- Return type:
- Warns:
UserWarning – Warns that this profile has no alternative parameterization.
Examples
>>> from socca.models import Beta >>> beta = Beta(xc=180.5, yc=45.2) >>> beta_refactored = beta.refactor()
- class socca.models.CustomProfile(parameters, profile, **kwargs)[source]¶
Bases:
ProfileUser-defined custom surface brightness profile.
Allows users to define arbitrary profile functions with custom parameters, enabling modeling of non-standard surface brightness distributions.
- Parameters:
parameters (list of dict) –
List of parameter specifications. Each dict should contain:
’name’: str, parameter name
’unit’: str, optional, physical unit (default: ‘not specified’)
’description’: str, optional, parameter description
profile (callable) – Function defining the profile. Should have signature profile(r, **params) where r is the elliptical radius and **params are the custom parameters.
**kwargs (dict) – Standard profile parameters (xc, yc, theta, e, cbox, positive).
Notes
The profile function is automatically JIT-compiled for performance.
All parameters in the parameters list are initialized to None and must be set before use.
The profile function should be compatible with JAX (use jax.numpy operations).
Examples
>>> from socca.models import CustomProfile >>> import jax.numpy as jp >>> # Define a custom Gaussian profile >>> def gaussian_profile(r, amplitude, sigma): ... return amplitude * jp.exp(-0.5 * (r / sigma)**2) >>> params = [ ... {'name': 'amplitude', 'unit': 'image', 'description': 'Peak value'}, ... {'name': 'sigma', 'unit': 'deg', 'description': 'Gaussian width'} ... ] >>> profile = CustomProfile(params, gaussian_profile, xc=180.5, yc=45.2) >>> profile.amplitude = 100.0 >>> profile.sigma = 0.01
- __init__(parameters, profile, **kwargs)[source]¶
Initialize a custom profile with user-defined parameters and function.
- Parameters:
parameters (list of dict) – Parameter specifications with ‘name’, ‘unit’, and ‘description’.
profile (callable) – Profile function with signature profile(r, **params).
**kwargs (dict) – Standard profile parameters (xc, yc, theta, e, cbox).
- class socca.models.Beta(**kwargs)[source]¶
Bases:
ProfileBeta profile for modeling galaxy clusters and elliptical galaxies.
The Beta profile describes a power-law density distribution commonly used for X-ray and radio observations of galaxy clusters. It has the form I(r) = Ic * (1 + (r/rc)^2)^(-beta).
- static profile(r, Ic, rc, alpha, beta)[source]¶
Beta profile surface brightness distribution.
The Beta profile describes a power-law density distribution commonly used for modeling galaxy clusters and elliptical galaxies.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Ic (float) – Central surface brightness (same units as image).
rc (float) – Core radius in degrees.
alpha (float) – Radial exponent parameter (default 2.0 for standard Beta profile).
beta (float) – Slope parameter (typically 0.4-1.0 for galaxy clusters).
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The generalized Beta profile is defined as: I(r) = Ic * [1 + (r/rc)^alpha]^(-beta)
With alpha=2, this reduces to the standard Beta profile commonly used in X-ray astronomy for modeling hot gas in galaxy clusters.
References
Cavaliere, A., & Fusco-Femiano, R. 1976, A&A, 49, 137
Examples
>>> import jax.numpy as jp >>> from socca.models import Beta >>> r = jp.linspace(0, 0.1, 100) >>> I = Beta.profile(r, Ic=100.0, rc=0.01, alpha=2.0, beta=0.5)
- class socca.models.gNFW(**kwargs)[source]¶
Bases:
ProfileGeneralized Navarro-Frenk-White (gNFW) profile.
The gNFW profile is a flexible three-parameter model commonly used to describe the surface brightness distribution of galaxy clusters. It generalizes the NFW profile with adjustable inner (gamma), intermediate (alpha), and outer (beta) slopes.
By default the projected surface brightness is computed via numerical Abel integration on a 1D radial grid that is then interpolated onto the image pixels (accurate and fast). Direct per-pixel integration is also available, as well as a pre-trained MLP emulator.
- Parameters:
emulator (bool or str or None, optional) – Controls whether to use the MLP emulator instead of numerical integration. Pass
Trueto use the bundled pre-trained model, a file-system path (strorPath) to load a custom checkpoint, orNone(default) to use numerical integration. Requiresflaxto be installed.interpolate (bool, optional) – Selects the evaluation strategy. If
True(default) the profile is computed on a 1D radial grid and the result is interpolated to the image pixels. IfFalsethe profile is evaluated directly at every pixel (slower but free from interpolation error). Applies to both the numerical integration and emulator paths.rz (array-like, optional) – Radial grid (in units of rc) for the interpolation path. Only used when
interpolateis True. Default is logspace(-7, 2, 1000).eps (float, optional) – Absolute and relative quadrature tolerance. Only used when
emulatoris None. Default is 1e-8.
- property alpha¶
- property beta¶
- property gamma¶
- class socca.models.Power(**kwargs)[source]¶
Bases:
ProfilePower law profile for surface brightness modeling.
The Power profile describes a simple power-law distribution of the form I(r) = Ic * (r/rc)^alpha.
- static profile(r, Ic, rc, alpha)[source]¶
Power law surface brightness distribution.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Ic (float) – Characteristic surface brightness (same units as image).
rc (float) – Scale radius in degrees.
alpha (float) – Power law slope parameter.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The Power law profile is defined as: I(r) = Ic * (r/rc)^(-alpha)
For positive alpha values, the profile decreases with radius.
Examples
>>> import jax.numpy as jp >>> from socca.models import Power >>> r = jp.linspace(0.001, 0.1, 100) >>> I = Power.profile(r, Ic=100.0, rc=0.01, alpha=2.0)
- class socca.models.TopHat(**kwargs)[source]¶
Bases:
ProfileTop-hat surface brightness profile.
The Top-hat profile describes a uniform surface brightness distribution within a cutoff radius: I(r) = 1 for
|r|< rc, and I(r) = 0 otherwise. This profile is useful for modeling flat-topped emission regions.- static profile(r, rc, Ic)[source]¶
Top-hat surface brightness distribution.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
rc (float) – Cutoff radius in degrees.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The Top-hat profile is defined as:
I(r) = 1 for
|r|< rc, and I(r) = 0 otherwiseThis profile produces a flat, uniform emission within the cutoff radius and zero emission outside.
- class socca.models.Sersic(**kwargs)[source]¶
Bases:
ProfileSersic profile for modeling elliptical galaxies and bulges.
The Sersic profile is a generalization of de Vaucouleurs’ law that describes the light distribution in elliptical galaxies and galactic bulges. The profile shape is controlled by the Sersic index (ns), with ns=1 corresponding to an exponential disk and ns=4 to a de Vaucouleurs profile.
- property ns¶
- static profile(r, Ie, re, ns)[source]¶
Sersic profile surface brightness distribution.
The Sersic profile is a generalization of de Vaucouleurs’ law and describes the light distribution in elliptical galaxies and bulges.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Ie (float) – Surface brightness at the effective radius (same units as image).
re (float) – Effective (half-light) radius in degrees.
ns (float) –
Sersic index (concentration parameter). Typical values:
ns = 0.5-1: disk-like profiles
ns = 4: de Vaucouleurs profile (elliptical galaxies)
ns > 4: highly concentrated profiles
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The Sersic profile is defined as: I(r) = Ie * exp(-bn * [(r/re)^(1/ns) - 1])
where bn is chosen such that re encloses half the total light. The parameter bn is approximated numerically and interpolated.
Common special cases:
ns = 1: Exponential profile
ns = 4: de Vaucouleurs profile (elliptical galaxies)
The valid range for ns is approximately 0.25 to 10.
References
Sersic, J. L. 1968, Atlas de Galaxias Australes Ciotti, L., & Bertin, G. 1999, A&A, 352, 447
Examples
>>> import jax.numpy as jp >>> from socca.models import Sersic >>> r = jp.linspace(0, 0.1, 100) >>> # de Vaucouleurs profile for elliptical galaxy >>> I = Sersic.profile(r, Ie=50.0, re=0.02, ns=4.0)
- class socca.models.Gaussian(**kwargs)[source]¶
Bases:
ProfileGaussian surface brightness profile.
The Gaussian profile describes a surface brightness distribution that follows a Gaussian function of radius: I(r) = Is * exp(-0.5 * (r/rs)^2). This profile is equivalent to a Sersic profile with index n=0.5.
- __init__(**kwargs)[source]¶
Initialize a Gaussian profile component.
- Parameters:
**kwargs (dict) –
Keyword arguments including:
- rsfloat, optional
Scale radius (standard deviation) in degrees.
- Isfloat, optional
Central surface brightness (same units as image).
- xc, ycfloat, optional
Centroid coordinates (inherited from Profile).
- thetafloat, optional
Position angle in radians (inherited from Profile).
- efloat, optional
Projected axis ratio (inherited from Profile).
- cboxfloat, optional
Projected boxiness (inherited from Profile).
- static profile(r, Is, rs)[source]¶
Gaussian surface brightness distribution.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Is (float) – Central surface brightness (same units as image).
rs (float) – Scale radius (standard deviation) in degrees.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The Gaussian profile is defined as:
I(r) = Is * exp(-0.5 * (r / rs)^2)
This is equivalent to a Sersic profile with n=0.5. The scale radius rs corresponds to the standard deviation of the Gaussian, and the half-width at half-maximum (HWHM) is approximately 1.177 * rs.
- class socca.models.Exponential(**kwargs)[source]¶
Bases:
ProfileExponential disk profile.
The exponential profile (Sersic index n=1) is the standard model for disk galaxies. It describes a surface brightness distribution that decays exponentially with radius: I(r) = Is * exp(-r/rs).
- static profile(r, Is, rs)[source]¶
Exponential disk profile surface brightness distribution.
The exponential profile (Sersic index n=1) describes the light distribution in disk galaxies and spiral arms.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Is (float) – Central surface brightness (same units as image).
rs (float) – Scale radius (scale length) in degrees.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The exponential profile is defined as: I(r) = Is * exp(-r / rs)
This is a special case of the Sersic profile with n=1 and is the canonical profile for modeling disk galaxies. The scale radius rs contains about 63% of the total light within that radius.
The exponential profile has no finite effective radius; the half-light radius is approximately 1.678 * rs.
References
Freeman, K. C. 1970, ApJ, 160, 811 van der Kruit, P. C., & Searle, L. 1981, A&A, 95, 105
Examples
>>> import jax.numpy as jp >>> from socca.models import Exponential >>> r = jp.linspace(0, 0.1, 100) >>> I = Exponential.profile(r, Is=100.0, rs=0.02)
- class socca.models.PolyExponential(**kwargs)[source]¶
Bases:
ExponentialExponential profile with polynomial modulation.
Extended exponential profile that includes polynomial terms to model deviations from a pure exponential disk. Based on the formalism from Mancera Pina et al., A&A, 689, A344 (2024).
- static profile(r, Is, rs, c1, c2, c3, c4, rc)[source]¶
Polynomial-exponential profile surface brightness distribution.
An exponential profile modulated by a 4th-order polynomial, providing flexibility to model deviations from pure exponential profiles in disk galaxies.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Is (float) – Central surface brightness (same units as image).
rs (float) – Exponential scale radius in degrees.
c1 (float) – Polynomial coefficients for 1st through 4th order terms.
c2 (float) – Polynomial coefficients for 1st through 4th order terms.
c3 (float) – Polynomial coefficients for 1st through 4th order terms.
c4 (float) – Polynomial coefficients for 1st through 4th order terms.
rc (float) – Reference radius for polynomial terms in degrees.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The profile is defined as:
- I(r) = Is * exp(-r/rs) * [1 + c1*(r/rc) + c2*(r/rc)^2
c3*(r/rc)^3 + c4*(r/rc)^4]
This profile allows modeling of:
Truncations or drops in outer regions
Central enhancements or deficits
Breaks in disk profiles
Type I, II, and III disk breaks
The polynomial modulation is normalized to the reference radius rc, which should typically be comparable to the scale radius rs.
References
Mancera Pina et al., A&A, 689, A344 (2024)
Examples
>>> import jax.numpy as jp >>> from socca.models import PolyExponential >>> r = jp.linspace(0, 0.1, 100) >>> # Profile with slight central enhancement >>> I = PolyExponential.profile(r, Is=100.0, rs=0.02, c1=0.1, ... c2=-0.05, c3=0.0, c4=0.0, rc=0.02)
- class socca.models.PolyExpoRefact(**kwargs)[source]¶
Bases:
ExponentialRefactored polynomial exponential profile with intensity coefficients.
Alternative parameterization of the polynomial exponential profile using intensity coefficients instead of polynomial coefficients. Based on the formalism from Mancera Pina et al., A&A, 689, A344 (2024).
- static profile(r, Is, rs, I1, I2, I3, I4, rc)[source]¶
Refactored polynomial-exponential profile using intensity coefficients.
Alternative parameterization of the polynomial-exponential profile using intensity coefficients rather than dimensionless polynomial coefficients.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Is (float) – Central surface brightness (same units as image).
rs (float) – Exponential scale radius in degrees.
I1 (float) – Intensity coefficients for polynomial terms (same units as Is).
I2 (float) – Intensity coefficients for polynomial terms (same units as Is).
I3 (float) – Intensity coefficients for polynomial terms (same units as Is).
I4 (float) – Intensity coefficients for polynomial terms (same units as Is).
rc (float) – Reference radius for polynomial terms in degrees.
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The profile is defined as:
- I(r) = exp(-r/rs) * [Is + I1*(r/rc) + I2*(r/rc)^2
I3*(r/rc)^3 + I4*(r/rc)^4]
This is mathematically equivalent to PolyExponential but uses intensity coefficients Ii instead of dimensionless coefficients ci. The relation is:
Ii = ci * Is
This parameterization may be more intuitive when the polynomial terms represent physical contributions with intensity units.
References
Mancera Pina et al., A&A, 689, A344 (2024)
See also
PolyExponentialAlternative parameterization with dimensionless coefficients.
- refactor()[source]¶
Convert to equivalent PolyExponential parameterization.
Transforms the intensity-based parameterization (I1, I2, I3, I4) to the dimensionless coefficient parameterization (c1, c2, c3, c4) of PolyExponential.
- Returns:
Equivalent profile with dimensionless coefficients ci = Ii / Is.
- Return type:
Notes
This conversion preserves the exact same surface brightness profile but expresses it using normalized polynomial coefficients.
Examples
>>> from socca.models import PolyExpoRefact >>> prof = PolyExpoRefact(Is=100, I1=10, I2=5, I3=0, I4=0) >>> prof_equiv = prof.refactor() >>> # prof_equiv is PolyExponential with c1=0.1, c2=0.05, c3=0, c4=0
- class socca.models.ModExponential(**kwargs)[source]¶
Bases:
ExponentialModified exponential profile with power-law modulation.
An exponential profile modulated by a power law to provide additional flexibility for modeling disk profiles with deviations from pure exponential behavior. Based on the formalism from Mancera Pina et al., A&A, 689, A344 (2024).
- static profile(r, Is, rs, rm, alpha)[source]¶
Generate modified exponential profile with power-law modulation.
An exponential profile modulated by a power law, providing an additional degree of freedom to model disk profiles with deviations from pure exponential behavior.
- Parameters:
r (ndarray) – Elliptical radius in degrees.
Is (float) – Central surface brightness (same units as image).
rs (float) – Exponential scale radius in degrees.
rm (float) – Modification radius where power law becomes important (deg).
alpha (float) – Power-law exponent (positive for enhancement, negative for suppression).
- Returns:
Surface brightness at radius r.
- Return type:
ndarray
Notes
The profile is defined as: I(r) = Is * exp(-r/rs) * (1 + r/rm)^alpha
This profile can model:
Truncations (alpha < 0)
Enhancements in outer regions (alpha > 0)
Smooth transitions between different slopes
The modification becomes significant at r ~ rm. For r << rm, the profile is approximately exponential. For r >> rm, the behavior depends on the sign of alpha.
References
Mancera Pina et al., A&A, 689, A344 (2024)
Examples
>>> import jax.numpy as jp >>> from socca.models import ModExponential >>> r = jp.linspace(0, 0.1, 100) >>> # Profile with outer truncation >>> I = ModExponential.profile(r, Is=100.0, rs=0.02, rm=0.05, alpha=-2.0)
- class socca.models.Point(**kwargs)[source]¶
Bases:
ComponentPoint source model for unresolved sources.
Models point sources (stars, quasars, unresolved AGN) that are handled efficiently in Fourier space. The source is convolved with the PSF to produce the observed image.
- Parameters:
**kwargs (dict) –
Keyword arguments including:
- xcfloat, optional
Right ascension in degrees.
- ycfloat, optional
Declination in degrees.
- Icfloat, optional
Integrated flux or peak brightness (same units as image).
- positivebool, optional
Whether to enforce positivity constraint.
- Variables:
xc (float) – Source right ascension (deg).
yc (float) – Source declination (deg).
Ic (float) – Source intensity.
Notes
Point sources are special-cased in the model evaluation:
Handled in Fourier space using phase shifts
Always PSF-convolved (not meaningful without PSF)
More efficient than evaluating a very narrow Gaussian
Can account for instrumental response at the source position
- __init__(**kwargs)[source]¶
Initialize a point source component.
- Parameters:
**kwargs (dict) – Keyword arguments including xc, yc, Ic, positive.
- static profile(xc, yc, Ic)[source]¶
Point source profile placeholder (not used).
Point sources are handled specially in Fourier space and do not use a standard radial profile function.
- Parameters:
xc (float) – Right ascension (deg).
yc (float) – Declination (deg).
Ic (float) – Source intensity.
Notes
This method is a placeholder to maintain API consistency. Point sources are actually evaluated in getmap() using Fourier space phase shifts.
- static getgrid()[source]¶
Point source grid placeholder (not used).
Point sources do not use a spatial grid; they are handled in Fourier space.
Notes
This method is a placeholder to maintain API consistency with other profile types that use getgrid() for coordinate transformations.
- getmap(img, convolve=False)[source]¶
Generate point source image via Fourier space phase shifts.
Creates a point source image by computing the appropriate phase shift in Fourier space and optionally convolving with the PSF.
- Parameters:
img (Image) – Image object containing FFT information and PSF.
convolve (bool, optional) – If True, convolve with PSF. If False, return unconvolved point source (delta function on pixel grid). Default is False.
- Returns:
Point source image on the image grid.
- Return type:
ndarray
- Raises:
ValueError – If any parameter is a prior distribution or set to None.
- Warns:
UserWarning – If convolve=True but no PSF is defined.
Notes
The point source is created using the Fourier shift theorem:
The source is placed at (xc, yc) via phase shifts
Multiplication by PSF in Fourier space performs convolution
More efficient than spatial convolution for point sources
The ‘pulse’ factor accounts for Fourier normalization
For point sources, PSF convolution is typically essential since an unconvolved point source is a delta function (single bright pixel).
Examples
>>> from socca.models import Point >>> from socca.data import Image >>> point = Point(xc=180.5, yc=45.2, Ic=1000.0) >>> img = Image('observation.fits') >>> psf_convolved = point.getmap(img, convolve=True)
- class socca.models.Background(**kwargs)[source]¶
Bases:
ComponentPolynomial background model for large-scale gradients.
Models smooth background variations using a 2D polynomial up to 3rd order. Useful for fitting sky background, instrumental gradients, or scattered light.
- Parameters:
**kwargs (dict) –
Keyword arguments including:
- rsfloat, optional
Reference radius for normalizing polynomial terms (deg).
- a0float, optional
Constant (0th order) term.
- a1x, a1yfloat, optional
Linear (1st order) terms in x and y.
- a2xx, a2xy, a2yyfloat, optional
Quadratic (2nd order) terms.
- a3xxx, a3xxy, a3xyy, a3yyyfloat, optional
Cubic (3rd order) terms.
- positivebool, optional
Whether to enforce positivity constraint.
- Variables:
rs (float) – Reference radius for polynomial normalization (deg).
a3yyy (a0, a1x, a1y, a2xx, a2xy, a2yy, a3xxx, a3xxy, a3xyy,) – Polynomial coefficients up to 3rd order.
Notes
The background is defined as:
- B(x,y) = a0 + a1x·x’ + a1y·y’
a2xx·x’^2 + a2xy·x’·y’ + a2yy·y’^2
a3xxx·x’^3 + a3xxy·x’^2·y’ + a3xyy·x’·y’^2 + a3yyy·y’^3
where x’ = x/rs and y’ = y/rs are normalized coordinates.
Coordinates are relative to the field center (CRVAL1, CRVAL2)
x coordinate includes cos(dec) correction for spherical geometry
Background is not PSF-convolved (assumed to vary on large scales)
Typically use low-order terms (0th and 1st) to avoid overfitting
- __init__(**kwargs)[source]¶
Initialize a polynomial background component.
- Parameters:
**kwargs (dict) – Keyword arguments including rs and polynomial coefficients a0, a1x, etc.
- static profile(x, y, a0, a1x, a1y, a2xx, a2xy, a2yy, a3xxx, a3xxy, a3xyy, a3yyy, rs)[source]¶
Evaluate 2D polynomial background on coordinate grids.
- Parameters:
x (ndarray) – Coordinate grids in degrees (relative to field center).
y (ndarray) – Coordinate grids in degrees (relative to field center).
a0 (float) – Constant term.
a1x (float) – Linear coefficients.
a1y (float) – Linear coefficients.
a2xx (float) – Quadratic coefficients.
a2xy (float) – Quadratic coefficients.
a2yy (float) – Quadratic coefficients.
a3xxx (float) – Cubic coefficients.
a3xxy (float) – Cubic coefficients.
a3xyy (float) – Cubic coefficients.
a3yyy (float) – Cubic coefficients.
rs (float) – Reference radius for normalization (deg).
- Returns:
Background values on the coordinate grid.
- Return type:
ndarray
Notes
Evaluates the polynomial:
- B = a0 + a1x·x’ + a1y·y’ + a2xx·x’^2 + a2xy·x’·y’ + a2yy·y’^2
a3xxx·x’^3 + a3xxy·x’^2·y’ + a3xyy·x’·y’^2 + a3yyy·y’^3
where x’ = x/rs and y’ = y/rs.
The normalization by rs keeps coefficients of different orders at comparable scales and improves numerical conditioning.
- static getgrid(grid, xc, yc)[source]¶
Compute relative coordinates for background evaluation.
Converts absolute celestial coordinates to coordinates relative to the field center, with spherical geometry correction.
- Parameters:
grid (Grid) – Grid object with .x and .y coordinate arrays (deg).
xc (float) – Reference RA (field center) in degrees.
yc (float) – Reference Dec (field center) in degrees.
- Returns:
xgrid (ndarray) – RA offsets in degrees (with cos(dec) correction).
ygrid (ndarray) – Dec offsets in degrees.
Notes
The RA offset includes a cos(dec) factor to account for spherical geometry, ensuring that distances are approximately correct on the sky.
- getmap(img, **kwargs)[source]¶
Generate background map on the image grid.
Evaluates the polynomial background across the entire image field.
- Parameters:
img (Image) – Image object containing grid and WCS information.
**kwargs (dict) – Ignored keyword arguments (for API compatibility).
- Returns:
Background map on the image grid.
- Return type:
ndarray
- Raises:
ValueError – If any parameter is a prior distribution or set to None.
- Warns:
UserWarning – If ‘convolve’ argument is provided (background is never convolved).
Notes
Background is evaluated relative to field center (CRVAL1, CRVAL2)
Background is not PSF-convolved (assumed to vary smoothly)
All polynomial coefficients must have fixed values
Examples
>>> from socca.models import Background >>> from socca.data import Image >>> bkg = Background(a0=10.0, a1x=0.1, a1y=-0.05, rs=1.0) >>> img = Image('observation.fits') >>> bkg_map = bkg.getmap(img)
- class socca.models.Disk(radial=<socca.models.radial.sersic.Sersic object>, vertical=<socca.models.disk.vertical.HyperSecantHeight object>, **kwargs)[source]¶
Bases:
Component3D disk model combining radial and vertical profiles with inclination.
Creates a realistic disk galaxy model by combining a radial surface brightness profile with a vertical density distribution. The model accounts for disk inclination via line-of-sight integration.
- Parameters:
radial (Profile, optional) – Radial profile defining in-plane surface brightness distribution. Default is Sersic(). Can be any 2D profile (Beta, Exponential, etc.).
vertical (Height, optional) – Vertical profile defining scale height and vertical structure. Default is HyperSecantHeight(). Can be HyperSecantHeight or ExponentialHeight.
**kwargs (dict) – Additional keyword arguments passed to Component.
- Variables:
Notes
The 3D disk model: 1. Defines density as rho(r,z) = radial_profile(r) * vertical_profile(z) 2. Integrates along the line of sight to obtain projected surface brightness 3. Accounts for inclination by rotating the disk before integration 4. Uses numerical integration (trapezoidal rule) along z-direction
The inclination angle is stored in vertical.inc:
inc = 0: face-on (looking down on the disk)
inc = pi/2: edge-on (viewing disk from the side)
Line-of-sight integration parameters (vertical.losdepth, vertical.losbins) control accuracy and computation time. Larger losdepth and more losbins improve accuracy for highly inclined disks.
Examples
>>> from socca.models import Disk, Exponential >>> from socca.models.disk.vertical import ExponentialHeight >>> # Create edge-on exponential disk >>> radial = Exponential(xc=180.5, yc=45.2, rs=0.02, Is=100) >>> vertical = ExponentialHeight(zs=0.002, inc=np.pi/2) >>> disk = Disk(radial=radial, vertical=vertical)
- __init__(radial=<socca.models.radial.sersic.Sersic object>, vertical=<socca.models.disk.vertical.HyperSecantHeight object>, **kwargs)[source]¶
Initialize a 3D disk model from radial and vertical profiles.
- Parameters:
Notes
The component ID is synchronized between the Disk, radial, and vertical components to ensure consistent parameter naming in composite models.
- getmap(img, convolve=False)[source]¶
Generate disk image via 3D line-of-sight integration.
Computes the projected surface brightness by integrating the 3D density distribution along the line of sight, accounting for inclination.
- Parameters:
img (Image) – Image object containing grid, PSF, and WCS information.
convolve (bool, optional) – If True, convolve the model with the PSF. Default is False.
- Returns:
Projected disk image on the image grid.
- Return type:
ndarray
- Raises:
ValueError – If any parameter is a prior distribution or set to None.
- Warns:
UserWarning – If convolve=True but no PSF is defined.
Notes
The algorithm: 1. Creates 3D coordinate grids (r, z) accounting for inclination 2. Evaluates 3D density rho(r, z) = radial(r) * vertical(z) 3. Integrates along line of sight using trapezoidal rule 4. Averages over sub-pixel sampling 5. Optionally convolves with PSF
Integration accuracy controlled by:
vertical.losdepth: extent of integration
vertical.losbins: number of integration points
For edge-on disks, use larger losdepth (≥ 5*zs) and more losbins (≥ 200).
Examples
>>> from socca.models import Disk, Exponential >>> from socca.models.disk.vertical import ExponentialHeight >>> import numpy as np >>> from socca.data import Image >>> radial = Exponential(xc=180.5, yc=45.2, rs=0.02, Is=100) >>> vertical = ExponentialHeight(zs=0.002, inc=np.pi/4) # 45 deg >>> disk = Disk(radial=radial, vertical=vertical) >>> img = Image('observation.fits') >>> disk_map = disk.getmap(img, convolve=True)
- static getgrid(grid, xc, yc, losdepth, losbins=200, theta=0.0, inc=0.0)[source]¶
Compute 3D disk coordinates with rotation and inclination.
Generates 3D coordinate grids (r, z) for disk model evaluation, accounting for position angle and inclination transformations needed for line-of-sight integration through an inclined disk.
- Parameters:
grid (Grid) – Grid object with .x and .y celestial coordinate arrays (deg).
xc (float) – Right ascension of disk center (deg).
yc (float) – Declination of disk center (deg).
losdepth (float) – Half-extent of line-of-sight integration (deg).
losbins (int, optional) – Number of integration points along line of sight. Default is 200.
theta (float, optional) – Position angle of disk major axis, east from north (rad). Default is 0.
inc (float, optional) – Inclination angle (0 = face-on, pi/2 = edge-on) (rad). Default is 0.
- Returns:
rcube (ndarray) – 4D array of radial distances from disk center (deg). Shape: (ssize, losbins, ysize, xsize).
zcube (ndarray) – 4D array of heights above/below disk midplane (deg). Shape: (ssize, losbins, ysize, xsize).
Notes
The transformation sequence: 1. Center coordinates on (xc, yc) 2. Apply spherical geometry correction (cos(dec)) 3. Rotate by position angle theta 4. Create line-of-sight grid from -losdepth to +losdepth 5. Apply inclination rotation 6. Compute cylindrical radius r and height z
The inclination is measured from face-on:
inc = 0: face-on, z-axis points toward observer
inc = pi/2: edge-on, disk in plane of sky
The 4D arrays allow vectorized evaluation of the 3D density function across the entire image with sub-pixel sampling and line-of-sight integration.
This function is JIT-compiled with static losbins for performance.
- parameters()[source]¶
Print formatted table of disk parameters from radial and vertical components.
Displays parameters from both the radial and vertical profile components, organized as ‘radial.parameter’ and ‘vertical.parameter’. Separates regular parameters from hyperparameters (integration settings).
Notes
Output format:
Model parameters¶
radial.xc [deg] : value | Right ascension of centroid radial.yc [deg] : value | Declination of centroid radial.rs [deg] : value | Scale radius vertical.zs [deg] : value | Scale height vertical.inc [rad] : value | Inclination angle (0=face-on)
Hyperparameters¶
vertical.losdepth [deg] : value | Half line-of-sight extent vertical.losbins [] : value | Number of integration points
Hyperparameters control the numerical integration accuracy but are not fitted parameters.
Examples
>>> from socca.models import Disk, Exponential >>> from socca.models.disk.vertical import ExponentialHeight >>> import numpy as np >>> radial = Exponential(xc=180.5, yc=45.2, rs=0.02, Is=100) >>> vertical = ExponentialHeight(zs=0.002, inc=np.pi/4) >>> disk = Disk(radial=radial, vertical=vertical) >>> disk.parameters() Model parameters ================ radial.xc [deg] : 1.8050E+02 | Right ascension of centroid radial.yc [deg] : 4.5200E+01 | Declination of centroid ...
- parlist()[source]¶
Return list of parameter names from radial and vertical components.
Collects all parameters from both the radial and vertical profile components, with names prefixed as ‘radial.’ and ‘vertical.’.
- Returns:
Combined list of parameter names from radial and vertical components. Names are prefixed with component type (e.g., ‘radial.xc’, ‘vertical.zs’).
- Return type:
list of str
Notes
This method is used internally by Model.addcomponent() when adding a Disk component to a composite model, ensuring all parameters from both sub-components are registered.
Examples
>>> from socca.models import Disk, Exponential >>> from socca.models.disk.vertical import ExponentialHeight >>> radial = Exponential(xc=180.5, yc=45.2, rs=0.02, Is=100) >>> vertical = ExponentialHeight(zs=0.002, inc=0.5) >>> disk = Disk(radial=radial, vertical=vertical) >>> disk.parlist() ['radial.xc', 'radial.yc', 'radial.theta', 'radial.e', 'radial.cbox', 'radial.rs', 'radial.Is', 'vertical.losdepth', 'vertical.losbins', 'vertical.inc', 'vertical.zs']
- class socca.models.SimpleBridge(radial=<socca.models.radial.beta.Beta object>, parallel=<socca.models.radial.tophat.TopHat object>, **kwargs)[source]¶
Bases:
BridgeSimple bridge model with multiplicative profile combination.
The SimpleBridge combines radial and parallel profiles multiplicatively, producing a surface brightness distribution of the form:
I(r, z) = Is * f_radial(r) * f_parallel(z)
This creates structures where the emission is the product of the two independent profile functions.
- Parameters:
radial (Profile, optional) – Profile describing emission perpendicular to the bridge axis. Default is Beta().
parallel (Profile, optional) – Profile describing emission along the bridge axis. Default is TopHat(), producing a uniform distribution along the bridge length.
**kwargs (dict) – Additional keyword arguments passed to Bridge base class.
See also
BridgeBase class with parameter descriptions.
MesaBridgeAlternative with mesa-like profile combination.
Examples
>>> from socca.models import SimpleBridge >>> bridge = SimpleBridge() >>> bridge.parameters()
Model parameters¶
xc [deg] : None | Right ascension of bridge centroid yc [deg] : None | Declination of bridge centroid rs [deg] : None | Scale radius of the bridge profiles Is [image] : None | Scale intensity of the bridge theta [rad] : 0.0000E+00 | Position angle of bridge major axis e [] : 5.0000E-01 | Projected ellipticity of the bridge profiles radial.alpha [] : 2.0000E+00 | Radial exponent radial.beta [] : 5.5000E-01 | Slope parameter
- class socca.models.MesaBridge(radial=None, parallel=None, **kwargs)[source]¶
Bases:
BridgeMesa bridge model with harmonic mean profile combination.
The MesaBridge combines radial and parallel profiles using a harmonic mean, producing a mesa-like (flat-topped) surface brightness distribution:
I(r, z) = Is / (1/f_radial(r) + 1/f_parallel(z))
This creates smooth transitions between the flat central region and the declining edges, resembling a mesa or table-top shape.
- Parameters:
radial (Profile, optional) – Profile for perpendicular direction. Default is Beta with alpha=8.0 and beta=1.0 for steep edges.
parallel (Profile, optional) – Profile for parallel direction. Default is Power with alpha=8.0 for steep drop-off.
**kwargs (dict) – Additional keyword arguments passed to Bridge base class.
Notes
The default parameters are chosen to produce a characteristic mesa-like shape with steep edges, suitable for modeling intracluster bridges with relatively uniform central emission.
See also
BridgeBase class with parameter descriptions.
SimpleBridgeAlternative with multiplicative profile combination.
Examples
>>> from socca.models import MesaBridge >>> bridge = MesaBridge() >>> bridge.parameters()
References
Hincks, A. D., et al., MRANS, 510, 3335 (2022) https://scixplorer.org/abs/2022MNRAS.510.3335H/abstract
Model parameters¶
xc [deg] : None | Right ascension of bridge centroid yc [deg] : None | Declination of bridge centroid rs [deg] : None | Scale radius of the bridge profiles Is [image] : None | Scale intensity of the bridge theta [rad] : 0.0000E+00 | Position angle of bridge major axis e [] : 5.0000E-01 | Axis ratio of the bridge profiles radial.alpha [] : 8.0000E+00 | Radial exponent radial.beta [] : 1.0000E+00 | Slope parameter parallel.alpha [] : 8.0000E+00 | Power law slope