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

The class for performing Schwarzschild modelling of density and kinematic data from the theoretical point of view (without "observational errors"). More...

#include <schwarzschild.h>

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

Public Member Functions

 CSchwModelQuadOpt (COrbitLibrary *_orbits, const vectorSchwData &_data, const CBasicQuadraticOptimizationSolver< smNumType > *_solver, const CBasicOrbitFilteringFnc *_orbitPenaltyFnc, double _lambda, double _penaltyConsLin, double _penaltyConsQuad, double sm_maxWeight)
 
virtual ~CSchwModelQuadOpt ()
 
virtual std::string performModelling ()
 perform the actual work, and return statistical information on success or error description on failure
 
virtual std::string getStatistics () const
 computes and returns statistical information. More...
 
virtual size_t numVariables () const
 total number of variables in the linear matrix (number of orbits plus slack variables, etc)
 
virtual size_t numConstraints () const
 total number of constraints in the linear matrix (in all data objects plus anything else)
 
virtual smNumType getLinearMatrixElem (size_t i1, size_t i2) const
 return the elements of linear constraint matrix to pass to the solver through a separately defined matrix interface object
 
void calcCoefDeviation (size_t coefIndex, double *deviation, size_t *numOrbits, size_t *numOrbitsNonzero) const
 compute deviations of a model constraint from its required value, normalized by coefNormFactor More...
 
double calcCoefComposition (size_t coefIndex, const CBasicOrbitFilteringFnc *weightFnc) const
 returns the value of the given constraint obtained by summing the contributions of all orbits weighted with their weight in the model and with the value of the provided weight function
 
bool isCorrect () const
 whether the orbit library is consistent with data objects
 
- Public Member Functions inherited from smile::CBasicSchwModel
 CBasicSchwModel (COrbitLibrary *_orbits, const vectorSchwData &_data)
 constructors of descendant classes would also take a pointer to some kind of solver as a parameter
 

Protected Attributes

const
CBasicQuadraticOptimizationSolver
< smNumType > * 
solver
 pointer to an existing instance of linear/quadratic solver
 
const CBasicOrbitFilteringFncorbitPenaltyFnc
 pointer to a function which assigns penalties to some kinds of orbits, may be NULL
 
const double lambda
 regularization term in quadratic optimization function, may be zero (then a linear optimization is performed)
 
const double penaltyConsLin
 linear penalty coefficient for constraint violation; if zero, then it means infinity (i.e. constraints must be satisfied exactly, but the problem may turn out to be infeasible)
 
const double penaltyConsQuad
 quadratic penalty coefficient for constraint violation
 
const double maxWeight
 if nonzero, puts upper limit to orbit weights
 
size_t numVar
 
size_t numConstr
 dimensions of linear matrix, set up in the constructor
 
size_t numConstrData
 total number of constraints in all data objects, without the additional tolerance constraints in the case penaltyCons<0
 
std::vector< smNumTyperhs
 required r.h.s. values of linear eq. system
 
bool correct
 flag set in the constructor to ensure that orbit library is consistent with Schwarzschild data
 
- Protected Attributes inherited from smile::CBasicSchwModel
COrbitLibraryorbits
 pointer to the orbit library
 
vectorSchwData data
 array of pointers to Schwarzschild data objects (which contain various constraints in the model)
 

Additional Inherited Members

- Public Types inherited from smile::CBasicSchwModel
enum  MODELTYPE { SM_LINEAR, SM_QUADRATIC }
 

Detailed Description

The class for performing Schwarzschild modelling of density and kinematic data from the theoretical point of view (without "observational errors").

The problem is solved by linear or quadratic optimization, in which the cost function takes the contribution to penalty from orbitPenaltyFnc if it is supplied (i.e. one may favour or disfavour some kinds of orbits, such as long-axis tubes or all regular orbits), and – in the case of quadratic programming – also the contribution from squared sum of orbit weights normalized to some priors, with a given proportionality constant lambda which has the meaning of regularization (smoothing). In other words, in the case without any penalty functions, it tries to achieve the distribution of orbit weights as close as possible to the priors (most often, uniform weights).

There are also two option to allow from departures from exact solution (i.e. not all constraints may be satisfied to machine precision). We denote deviation of a constraint as $ dev[c] = abs| \sum_o M[o][c] sol[o] - rhs[c] | $.

In the first approach, when penaltyCons>0, a linear term is added into the cost function which is just $ \sum_c penaltyCons * dev[c] $.

In the second approach, when penaltyCons<0, another set of N_c equations is introduced to keep the deviation within a given tolerance interval (as in van den Bosch et al.2008): $ |dev[c]| < |penaltyCons * rhs[c]| $.

The entire cost function is written as

\[ \mathcal{L} = \lambda/N_o * \sum_o (sol[o]/prior[o])^2 + \sum_o penaltyOrb[o] * sol[o]/M_{total} + \sum_c penaltyCons[c] * dev[c] . \]

If the problem is solved by linear programming, the first terms is omitted. If penaltyCons<0, the last term is absent and the tolerance intervals are absorbed by the additional set of non-negative variables.

Presently, some simplifications are in effect: priors for orbit weights are taken to be uniform prior[o] = M_total/N_o, and penalty terms for constraint violation are also the same for all constraints, but normalized to the intrinsic "normalization factors" provided by corresponding data objects.

Constructor & Destructor Documentation

smile::CSchwModelQuadOpt::CSchwModelQuadOpt ( COrbitLibrary _orbits,
const vectorSchwData _data,
const CBasicQuadraticOptimizationSolver< smNumType > *  _solver,
const CBasicOrbitFilteringFnc _orbitPenaltyFnc,
double  _lambda,
double  _penaltyConsLin,
double  _penaltyConsQuad,
double  sm_maxWeight 
)
Parameters
_orbitspointer to an existing orbit library
_datapointers to the array of Schwarzschild data objects (density+kinematic data)
_solverpointer to an instance of quad.opt. solver. Will be deleted when the object is destroyed! If NULL then no modelling could be performed, but statistics may be computed anyway
_orbitPenaltyFncpointer to an instance of a function evaluating penalties for orbits, or NULL for no such penalty. Will be deleted when the object is destroyed!
_lambdacoefficients in the penalty function for quadratic optimization; if 0 then a linear optimization is performed
_penaltyConsLinpenalties for constraint violation; if 0 then constraints must be satisfied exactly or the model is infeasible
_penaltyConsQuadquadratic penalty term for constraint violation
sm_maxWeightupper bound on orbit weight
virtual smile::CSchwModelQuadOpt::~CSchwModelQuadOpt ( )
inlinevirtual

destroys solver and penaltyfnc, but doesn't touch orbits or schwData

Member Function Documentation

void smile::CSchwModelQuadOpt::calcCoefDeviation ( size_t  coefIndex,
double *  deviation,
size_t *  numOrbits,
size_t *  numOrbitsNonzero 
) const

compute deviations of a model constraint from its required value, normalized by coefNormFactor

Parameters
[in]coefIndexis the index of the constraint
[out]deviationreturns the difference between the actual and required value of this constraint (if it is small enough, then it is rounded to zero)
[out]numOrbitsreturns the number of orbits that contribute anything to the value of this constraint
[out]numOrbitsNonzeroreturns the number of orbits with non-zero weight that contribute to this constraint
std::string smile::CSchwModelQuadOpt::getStatistics ( ) const
virtual

computes and returns statistical information.

The information returned may describe, e.g. how much the model differs from data, etc. Uses existing Schwarzschild data and orbit weights, without re-doing modelling.

Implements smile::CBasicSchwModel.


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