![]() |
SMILE
v2.5
Schwarzschild Modelling Interactive expLoratory Environment
|
Base class defining density model (without corresponding potential). More...
#include <potential.h>

Public Types | |
| enum | POTENTIALTYPE { PT_UNKNOWN, PT_DIRECT, PT_COMPOSITE, PT_COEFS, PT_NB, PT_BSE, PT_BSECOMPACT, PT_SPLINE, PT_CYLSPLINE, PT_LOG, PT_HARMONIC, PT_SCALEFREE, PT_SCALEFREESH, PT_SPHERICAL, PT_DEHNEN, PT_MIYAMOTONAGAI, PT_FERRERS, PT_PLUMMER, PT_ISOCHRONE, PT_PERFECTELLIPSOID, PT_NFW, PT_SERSIC, PT_EXPDISK, PT_ELLIPSOIDAL, PT_MGE } |
| list of all existing types of density or density/potential models, each of them implemented in its own class More... | |
| enum | SYMMETRYTYPE { ST_NONE = 0, ST_REFLECTION = 1, ST_PLANESYM = 2, ST_ZROTSYM = 4, ST_SPHSYM = 8, ST_TRIAXIAL = ST_REFLECTION | ST_PLANESYM, ST_AXISYMMETRIC = ST_TRIAXIAL | ST_ZROTSYM, ST_SPHERICAL = ST_AXISYMMETRIC | ST_SPHSYM, ST_DEFAULT = ST_TRIAXIAL } |
| Type of symmetry. More... | |
Public Member Functions | |
| virtual POTENTIALTYPE | PotentialType () const =0 |
| enumerable potential type | |
| virtual const char * | PotentialName () const =0 |
| string representation of potential type | |
| virtual SYMMETRYTYPE | symmetry () const =0 |
| returns symmetry type of this potential | |
| virtual double | Rho (double X, double Y, double Z, double t=0) const =0 |
| returns density at given coordinates, this should obviously be overriden in derivative classes | |
| virtual double | Mass (const double r) const |
| returns mass inside given radius (approximately! not necessary to integrate density over sphere, just a rough estimate used e.g. in choosing radial nodes of Schwarzschild grid) | |
| virtual double | totalMass () const |
| returns estimated M(r=infinity) or -1 if mass is infinite | |
| double | getRadiusByMass (const double m) const |
| solves for Mass(r)=m | |
| void | getRadiiByMass (const vectord &masses, vectord *radii) const |
| solves for Mass(r)=m for an array of sorted values of m (more efficient than doing it one-by-one) | |
| bool | checkMassMonotonic () const |
| safety measure: check (roughly) that mass is increasing with radius | |
| bool | checkDensityNonzero () const |
| another safety measure: check that density doesn't drop to zero along any of three axes (important to assess spherical-harmonic approximation quality) | |
| virtual double | getGamma () const |
| returns inner density slope estimate (only used in BSE potential expansion for the automatic selection of shape parameter Alpha) | |
Private Member Functions | |
| CDensity & | operator= (const CDensity &) |
| assignment forbidden, use clone() for CPotential | |
Base class defining density model (without corresponding potential).
May be used for initializing density for general-purpose potential approximations (BSE and Spline).
list of all existing types of density or density/potential models, each of them implemented in its own class
| Enumerator | |
|---|---|
| PT_UNKNOWN |
undefined/default value |
| PT_DIRECT |
direct evaluation of potential from Poisson equaiton: CPotentialDirect |
| PT_COMPOSITE |
a superposition of multiple potential instances: CPotentialComposite |
| PT_COEFS |
not an actual density model, but a way to load pre-computed coefficients of potential expansion |
| PT_NB |
a set of frozen particles: CPotentialNB |
| PT_BSE |
basis-set expansion for infinite systems: CPotentialBSE |
| PT_BSECOMPACT |
basis-set expansion for systems with non-singular density and finite extent: CPotentialBSECompact |
| PT_SPLINE |
spline spherical-harmonic expansion: CPotentialSpline |
| PT_CYLSPLINE |
expansion in azimuthal angle with two-dimensional meridional-plane interpolating splines: CPotentialCylSpline |
| PT_LOG |
logaritmic potential: CPotentialLog |
| PT_HARMONIC |
simple harmonic oscillator: CPotentialHarmonic |
| PT_SCALEFREE |
single power-law density profile: CPotentialScaleFree |
| PT_SCALEFREESH |
spherical-harmonic approximation to a power-law density: CPotentialScaleFreeSH |
| PT_SPHERICAL |
arbitrary spherical mass model: CPotentialSpherical |
| PT_DEHNEN |
Dehnen(1993) density model: CPotentialDehnen. |
| PT_MIYAMOTONAGAI |
Miyamoto-Nagai(1975) flattened model: CPotentialMiyamotoNagai. |
| PT_FERRERS |
Ferrers finite-extent profile: CPotentialFerrers. |
| PT_PLUMMER |
Plummer model: CDensityPlummer. |
| PT_ISOCHRONE |
isochrone model: CDensityIsochrone |
| PT_PERFECTELLIPSOID |
Kuzmin/de Zeeuw integrable potential: CDensityPerfectEllipsoid. |
| PT_NFW |
Navarro-Frenk-White profile: CDensityNFW. |
| PT_SERSIC |
Sersic density profile: CDensitySersic. |
| PT_EXPDISK |
exponential (in R) disk with a choice of vertical density profile: CDensityExpDisk |
| PT_ELLIPSOIDAL |
a generalization of spherical mass profile with arbitrary axis ratios: CDensityEllipsoidal |
| PT_MGE |
Multi-Gaussian expansion: CDensityMGE. |
Type of symmetry.
used mainly in spherical-harmonic expansion potentials to compute and use only a subset of coefficients based on the symmetry properties of the underlying potential)
1.8.8