![]()  | 
  
    SMILE
    v2.5
    
   Schwarzschild Modelling Interactive expLoratory Environment 
   | 
 
angular expansion of potential in azimuthal angle with coefficients being 2d spline functions of R,z More...
#include <potential.h>


Public Member Functions | |
| CPotentialCylSpline (size_t _Ncoefs_R, size_t _Ncoefs_z, size_t _Ncoefs_phi, const CDensity *density, double radius_min=0, double radius_max=0, double z_min=0, double z_max=0) | |
| init potential from analytic mass model specified by its density profile (using CPotentialDirect for intermediate potential computation)  | |
| CPotentialCylSpline (size_t _Ncoefs_R, size_t _Ncoefs_z, size_t _Ncoefs_phi, const CPotential *potential, double radius_min=0, double radius_max=0, double z_min=0, double z_max=0) | |
| init potential from analytic mass model specified by a potential-density pair  | |
| CPotentialCylSpline (const vectord &gridR, const vectord &gridz, const std::vector< vectord > &coefs) | |
| init potential from stored coefficients  | |
| CPotentialCylSpline (size_t _Ncoefs_R, size_t _Ncoefs_z, size_t _Ncoefs_phi, const CPointMassSet< double > &points, SYMMETRYTYPE _sym=ST_TRIAXIAL, double radius_min=0, double radius_max=0, double z_min=0, double z_max=0) | |
| init potential from N-body snapshot  | |
| virtual CPotential * | clone () const | 
| Return a pointer to a copy of this instance of potential.  More... | |
| virtual double | Mass (const double r) const | 
| a modified mass estimator that only integrates withing the grid definition radius  | |
| virtual POTENTIALTYPE | PotentialType () const | 
| enumerable potential type  | |
| virtual const char * | PotentialName () const | 
| string representation of potential type  | |
| 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 | 
| Compute forces and, optionally, force derivatives at a given point.  More... | |
| void | getCoefs (vectord *gridR, vectord *gridz, std::vector< vectord > *coefs) const | 
| retrieve coefficients of potential approximation.  More... | |
  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 | |
| void | initPot (size_t _Ncoefs_R, size_t _Ncoefs_z, size_t _Ncoefs_phi, const CPotential *potential, double radius_min, double radius_max, double z_min, double z_max) | 
| create interpolation grid and compute potential at grid nodes  | |
| void | initSplines (const std::vector< vectord > &coefs) | 
| create 2d interpolation splines in scaled R,z grid  | |
| void | initDefault () | 
| called as a default initialization when everything else fails  | |
| double | computePhi_m (double R, double z, int m, const CPotential *potential) const | 
| compute m-th azimuthal harmonic of potential (either by Fourier transform or calling the corresponding method of CPotentialDirect)  | |
Private Attributes | |
| SYMMETRYTYPE | mysymmetry | 
| may have different type of symmetry  | |
| vectord | grid_R | 
| vectord | grid_z | 
| nodes of the grid in cylindrical radius and vertical direction  | |
| double | Rscale | 
| scaling coefficient for transforming the interpolated potential; computed as -Phi(0)/Mtotal.  | |
| std::vector< interp2d_spline * > | splines | 
| array of 2d splines (for each m-component in the expansion in azimuthal angle)  | |
| double | C00 | 
| double | C20 | 
| double | C40 | 
| double | C22 | 
| multipole coefficients for extrapolation beyond the grid  | |
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... | |
angular expansion of potential in azimuthal angle with coefficients being 2d spline functions of R,z
      
  | 
  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.
      
  | 
  virtual | 
Compute forces and, optionally, force derivatives at a given point.
| [in] | xyz | - coordinates of the point to compute forces (array of 3 numbers) | 
| [in] | t | - time to compute forces (matters only if potential is time-dependent) | 
| [out] | force | - computed values of -d Phi/d x_i; output array must exist and contain N_DIM values. | 
| [out] | forceDeriv | - if not NULL, then also compute the second derivatives of potential (which is a symmetric matrix of size N_DIM^2, thus the array must be of size N_DIM*(N_DIM+1)/2 ): first 3 values contain  , second three contain mixed derivatives  .  | 
Implements smile::CPotential.
| void smile::CPotentialCylSpline::getCoefs | ( | vectord * | gridR, | 
| vectord * | gridz, | ||
| std::vector< vectord > * | coefs | ||
| ) | const | 
retrieve coefficients of potential approximation.
| [out] | gridR | will be filled with the array of R-values of grid nodes | 
| [out] | gridz | will be filled with the array of z-values of grid nodes | 
| [out] | coefs | will contain array of sequentially stored 2d arrays (the size of the outer array equals the number of terms in azimuthal expansion, inner arrays contain gridR.size()*gridz.size() values). All output array parameters must point to existing arrays. | 
      
  | 
  virtual | 
Return potential at a given spatial point (possibly a time-varying one).
Implements smile::CPotential.
 1.8.8