SMILE  v2.5
Schwarzschild Modelling Interactive expLoratory Environment
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
smile::CPotentialSpline Class Reference

spherical-harmonic expansion of potential with coefficients being spline functions of radius More...

#include <potential.h>

Inheritance diagram for smile::CPotentialSpline:
Inheritance graph
[legend]
Collaboration diagram for smile::CPotentialSpline:
Collaboration graph
[legend]

Public Member Functions

template<typename NumT >
 CPotentialSpline (size_t _Ncoefs_radial, size_t _Ncoefs_angular, const CPointMassSet< NumT > &points, SYMMETRYTYPE _sym=ST_TRIAXIAL, double smoothfactor=0, const vectord *_gridradii=NULL)
 init potential from N point masses with given assumed symmetry type, may also provide desired grid radii (otherwise assigned automatically)
 
 CPotentialSpline (const vectord &_gridradii, const std::vector< vectord > &_coefs)
 init potential from stored SHE coefficients at given radii
 
 CPotentialSpline (size_t _Ncoefs_radial, size_t _Ncoefs_angular, const CDensity *density, const vectord *_gridradii=NULL)
 init potential from analytic mass model, may also provide desired grid radii (if not given, assign automatically)
 
template<typename NumT >
 CPotentialSpline (size_t _Ncoefs_angular, const CPointMassSet< NumT > &points, const vectord *radii, std::vector< vectord > *coefsArray)
 not really a constructor, but a way to compute potential expansion coefficients for a set of discrete point masses at given radii. More...
 
virtual CPotentialclone () const
 Return a pointer to a copy of this instance of potential. More...
 
virtual POTENTIALTYPE PotentialType () const
 enumerable potential type
 
virtual const char * PotentialName () const
 string representation of potential type
 
virtual double Mass (const double r) const
 faster estimate of M(r) using l=0 terms only
 
size_t getNcoefs_radial () const
 return the number of radial points in the spline (excluding r=0)
 
void getCoefs (vectord *radii, std::vector< vectord > *coefsArray, bool useNodes=true) const
 compute SHE coefficients at given radii. More...
 
- Public Member Functions inherited from smile::CPotentialSH
 CPotentialSH (size_t _Ncoefs_angular)
 
virtual SYMMETRYTYPE symmetry () const
 returns symmetry type of this potential
 
virtual double Rho (double X, double Y, double Z, double t=0) const
 returns density at given coordinates, this should obviously be overriden in derivative classes
 
virtual double Phi (double X, double Y, double Z, double t=0) const
 Return potential at a given spatial point (possibly a time-varying one). More...
 
virtual void Force (const double xyz[N_DIM], const double t, double *force, double *forceDeriv=NULL) const
 common function for all derivative classes
 
size_t getNcoefs_angular () const
 return l_max – the order of angular expansion
 
- Public Member Functions inherited from smile::CDensity
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)
 

Static Public Member Functions

static const char * myName ()
 

Private Member Functions

virtual void computeSHCoefs (const double r, double coefsF[], double coefsdFdr[], double coefsd2Fdr2[]) const
 reimplemented function to compute the coefficients of angular expansion of potential at the given radius
 
void checkSymmetry (const std::vector< vectord > &coefsArray)
 assigns symmetry class if some coefficients are (near-)zero. More...
 
void initDefault ()
 called as a default initialization when everything else fails
 
void initSpline (const vectord &radii, const std::vector< vectord > &coefsArray)
 create spline objects for all non-zero spherical harmonics from the supplied radii and coefficients. More...
 
template<typename NumT >
void computeCoefsFromPoints (const CPointMassSet< NumT > &points, const vectord *srcradii, vectord *outradii, std::vector< vectord > *outcoefs)
 calculate (non-smoothed) spherical harmonic coefs for a discrete point mass set. More...
 
template<typename NumT >
void prepareCoefsDiscrete (const CPointMassSet< NumT > &points, double smoothfactor, const vectord *userradii)
 create smoothing splines from the coefficients computed at each particle's radius. More...
 
void prepareCoefsAnalytic (const CDensity *density, const vectord *srcradii)
 compute expansion coefficients from an analytical mass profile. More...
 
void coef0 (double r, double *val, double *der, double *der2) const
 evaluate value and optionally up to two derivatives of l=0 coefficient, taking into account extrapolation beyond the grid definition range and log-scaling of splines
 
void coeflm (size_t lm, double r, double xi, double *val, double *der, double *der2, double c0val, double c0der=0, double c0der2=0) const
 evaluate value, and optionally first and second derivative of l>0 coefficients (lm is the combined index of angular harmonic >0); corresponding values for 0th coef must be known
 

Private Attributes

size_t Ncoefs_radial
 number of radial coefficients (excluding the one at r=0)
 
vectord gridradii
 defines nodes of radial grid in splines
 
double minr
 
double maxr
 definition range of splines; extrapolation beyond this radius
 
double gammain
 
double coefin
 slope and coef. for extrapolating potential inside minr (spherically-symmetric part, l=0)
 
double gammaout
 
double coefout
 
double der2out
 slope and coef. for extrapolating potential outside maxr (spherically-symmetric part, l=0)
 
double potcenter
 
double potmax
 
double potminr
 (abs.value) potential in the center (for transformation of l=0 spline), at the outermost spline node, and at 1st spline node
 
std::vector< gsl_spline * > splines
 spline coefficients at each harmonic
 
vectord slopein
 
vectord slopeout
 slope of coefs for l>0 for extrapolating inside rmin/outside rmax
 

Additional Inherited Members

- Public Types inherited from smile::CDensity
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...
 
- Protected Member Functions inherited from smile::CPotentialSH
void assignlmrange ()
 assigns the above variables based on mysymmetry, should be called whenever mysymmetry has changed
 
- Protected Attributes inherited from smile::CPotentialSH
size_t Ncoefs_angular
 l_max, the order of angular expansion (0 means spherically symmetric model)
 
SYMMETRYTYPE mysymmetry
 may have different type of symmetry
 
int lmax
 
int lstep
 
int mmin
 
int mmax
 
int mstep
 range of angular coefficients used for given symmetry
 

Detailed Description

spherical-harmonic expansion of potential with coefficients being spline functions of radius

Constructor & Destructor Documentation

template<typename NumT >
smile::CPotentialSpline::CPotentialSpline ( size_t  _Ncoefs_angular,
const CPointMassSet< NumT > &  points,
const vectord radii,
std::vector< vectord > *  coefsArray 
)

not really a constructor, but a way to compute potential expansion coefficients for a set of discrete point masses at given radii.

Parameters
[in]pointscontains the source set of point masses;
[in]_Ncoefs_angularis the order of angular expansion (l_max);
[in]radiicontains the array of radii at which the coefficients should be computed;
[out]coefsArraycontains the computed coefficients. Smoothing is not performed, splines are not initialized and the potential is not usable after this call!

Member Function Documentation

void smile::CPotentialSpline::checkSymmetry ( const std::vector< vectord > &  coefsArray)
private

assigns symmetry class if some coefficients are (near-)zero.

called on an intermediate step after computing coefs but before initializing splines

CPotential * smile::CPotentialSpline::clone ( ) const
virtual

Return a pointer to a copy of this instance of potential.

A standard copy constructor or assignment is disabled because of different amount of data needed to be copied in different derived classes).

Implements smile::CPotential.

template<typename NumT >
void smile::CPotentialSpline::computeCoefsFromPoints ( const CPointMassSet< NumT > &  points,
const vectord srcradii,
vectord outradii,
std::vector< vectord > *  outcoefs 
)
private

calculate (non-smoothed) spherical harmonic coefs for a discrete point mass set.

Parameters
[in]pointscontains the positions and masses of particles that serve as the input. There are two modes of operation, depending on the content of "outradii" vector:
  • in the context of initialization of potential from a discrete point set, "srcradii"=NULL; "outradii" will be initialized with radii of all particles sorted in ascending order, and "outcoefs" will contain coefficients of SH expansion for all points, "outradii" and "outcoefs" should point to existing arrays; in case of error outcoefs is set to zero length.
    Note
    : the order of indexes in outcoefs is swapped w.r.t. standard convention in this file: outcoefs[coefind][pointind] - this is done to save memory by not allocating arrays for unused coefind'ices.
  • in the context of SHgrid Schwarzschild modelling, this is called from a quasi-constructor for the purpose of computing SH coefs at a given set of radii. This set is passed as "srcradii" array, and coefs are returned in "outcoefs" in the standard way: outcoefs[pointind][coefind]. "outcoefs" should point to an existing array, in case of error it is filled with zeroes.
Parameters
[in]srcradiimay contain user-defined set of radii to compute the expansion coefficients.
[out]outradiipoints to an existing array to be filled with radial grid points.
[out]outcoefspoints to an existing 2d array that will contain computed coefficients.
void smile::CPotentialSpline::getCoefs ( vectord radii,
std::vector< vectord > *  coefsArray,
bool  useNodes = true 
) const

compute SHE coefficients at given radii.

Parameters
[in]useNodestells whether the radii in question are the original spline nodes (true) or user-specified radii (true);
[in,out]radiiif useNodes=true, this array is filled with spline nodes [output parameter, but must point to an existing array], otherwise the coefficients are computed at the supplied radii [input parameter];
[out]coefsArrayis filled by the computed coefficients (must point to an existing array)
void smile::CPotentialSpline::initSpline ( const vectord radii,
const std::vector< vectord > &  coefsArray 
)
private

create spline objects for all non-zero spherical harmonics from the supplied radii and coefficients.

radii should have "Ncoefs_radial" elements and coefsArray - "Ncoefs_radial * (Ncoefs_angular+1)^2" elements

void smile::CPotentialSpline::prepareCoefsAnalytic ( const CDensity density,
const vectord srcradii 
)
private

compute expansion coefficients from an analytical mass profile.

Parameters
[in]densityis the density profile to be approximated
[in]srcradiimay contain the radial grid; if NULL then it is assigned automatically.
template<typename NumT >
void smile::CPotentialSpline::prepareCoefsDiscrete ( const CPointMassSet< NumT > &  points,
double  smoothfactor,
const vectord userradii 
)
private

create smoothing splines from the coefficients computed at each particle's radius.

Parameters
[in]pointsis the array of particle coordinates and masses
[in]userradiimay contain the radial grid to construct the splines; if it is NULL then the radial grid is constructed automatically taking into account the range of radii covered by source points.
[in]smoothfactordetermines how much smoothing is applied to the spline (good results are obtained for smoothfactor=1-2).

!!FIXME: what happens if the outermost pointRadius is < outerBinRadius ?


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