SlowDown Algorithmic Regularization (SDAR)
Algorithmic Regularization with slowdown method for integrating few-body motions
AR::Information< Tparticle, Tpcm > Class Template Reference

A class contains information (e.g. parameters, binary tree, indices) about the particle group. More...

#include <information.h>

Collaboration diagram for AR::Information< Tparticle, Tpcm >:
[legend]

Public Member Functions

 Information ()
 

a list of binary tree that contain the hierarchical orbital parameters of the particle group.

More...
 
bool checkParams ()
 check whether parameters values are correct initialized More...
 
void reserveMem (const int _nmax)
 reserve memory of binarytree list More...
 
BinaryTree< Tparticle > & getBinaryTreeRoot () const
 get the root of binary tree More...
 
Float calcDsElliptic (BinaryTree< Tparticle > &_bin, const Float &_G)
 
Float calcDsHyperbolic (BinaryTree< Tparticle > &_bin, const Float &_G)
 
void calcDsMinKeplerIter (Float &_ds_over_ebin_min_bin, Float &_ds_min_bin, Float &_ds_min_hyp, Float &_etot_sd, const Float &_G, const Float &_nest_sd_up, BinaryTree< Tparticle > &_bin, const int _intergrator_order)
 iteration function to calculate average kepler ds for a binary tree More...
 
Float calcDsKeplerBinaryTree (BinaryTree< Tparticle > &_bin, const int _int_order, const Float &_G, const Float &_ds_scale)
 calculate average kepler ds iterately for a binary tree More...
 
void calcDsAndStepOption (const int _int_order, const Float &_G, const Float &_ds_scale)
 calculate ds from the inner most binary with minimum period, determine the fix step option More...
 
void generateBinaryTree (COMM::ParticleGroup< Tparticle, Tpcm > &_particles, const Float _G)
 generate binary tree for the particle group More...
 
bool checkAndSetBinaryPairIDIter (BinaryTree< Tparticle > &_bin, const bool _reset_flag)
 check binary tree member pair id, if consisent, return ture. otherwise set the member pair id More...
 
void clear ()
 clear function More...
 
void printColumnTitle (std::ostream &_fout, const int _width=20)
 print titles of class members using column style More...
 
void printColumn (std::ostream &_fout, const int _width=20)
 print data of class members using column style More...
 
void writeBinary (FILE *_fout) const
 write class data to file with binary format More...
 
void readBinary (FILE *_fin)
 read class data to file with binary format More...
 

Public Attributes

Float ds
 
Float time_offset
 

initial step size for integration

More...
 
Float r_break_crit
 

offset of time to obtain real physical time (real time = TimeTransformedSymplecticIntegrator:time_ + info.time_offset)

More...
 
FixStepOption fix_step_option
 
COMM::List< BinaryTree< Tparticle > > binarytree
 

fix step option for integration

More...
 

Detailed Description

template<class Tparticle, class Tpcm>
class AR::Information< Tparticle, Tpcm >

A class contains information (e.g. parameters, binary tree, indices) about the particle group.

The member of this class should not be the data that must be recored and should be possible calculated any time based on the main class (TimeTransformedSymplecticIntegrator) data. This class must be inherited when a different Information class is applied in the main class, since the binary tree is used in the integration for slowdown factor The basic members are used in the integration are \ ds: integration step size \ binarytree: the Kepler orbital parameters of the hierarchical systems and slowdown factors \ fix_step_option: option to control whether the adjustment of step sizes are used \

Constructor & Destructor Documentation

◆ Information()

template<class Tparticle , class Tpcm >
AR::Information< Tparticle, Tpcm >::Information ( )
inline

a list of binary tree that contain the hierarchical orbital parameters of the particle group.

initializer, set ds to zero, fix_step_option to none

Member Function Documentation

◆ calcDsAndStepOption()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::calcDsAndStepOption ( const int  _int_order,
const Float _G,
const Float _ds_scale 
)
inline

calculate ds from the inner most binary with minimum period, determine the fix step option

Estimate ds first from the inner most binary orbit (eccentric anomaly), set fix_step_option to later

Parameters
[in]_int_orderaccuracy order of the symplectic integrator.
[in]_Ggravitational constant
[in]_ds_scalescaling factor to determine ds

◆ calcDsElliptic()

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::calcDsElliptic ( BinaryTree< Tparticle > &  _bin,
const Float _G 
)
inline
Here is the caller graph for this function:

◆ calcDsHyperbolic()

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::calcDsHyperbolic ( BinaryTree< Tparticle > &  _bin,
const Float _G 
)
inline
Here is the caller graph for this function:

◆ calcDsKeplerBinaryTree()

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::calcDsKeplerBinaryTree ( BinaryTree< Tparticle > &  _bin,
const int  _int_order,
const Float _G,
const Float _ds_scale 
)
inline

calculate average kepler ds iterately for a binary tree

use calcDsMinKeplerIter

Here is the caller graph for this function:

◆ calcDsMinKeplerIter()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::calcDsMinKeplerIter ( Float _ds_over_ebin_min_bin,
Float _ds_min_bin,
Float _ds_min_hyp,
Float _etot_sd,
const Float _G,
const Float _nest_sd_up,
BinaryTree< Tparticle > &  _bin,
const int  _intergrator_order 
)
inline

iteration function to calculate average kepler ds for a binary tree

Here is the caller graph for this function:

◆ checkAndSetBinaryPairIDIter()

template<class Tparticle , class Tpcm >
bool AR::Information< Tparticle, Tpcm >::checkAndSetBinaryPairIDIter ( BinaryTree< Tparticle > &  _bin,
const bool  _reset_flag 
)
inline

check binary tree member pair id, if consisent, return ture. otherwise set the member pair id

Parameters
[in]_binbinary tree to check
[in]_reset_flagif true, reset pair id to zero
Here is the caller graph for this function:

◆ checkParams()

template<class Tparticle , class Tpcm >
bool AR::Information< Tparticle, Tpcm >::checkParams ( )
inline

check whether parameters values are correct initialized

Returns
true: all correct

◆ clear()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::clear ( )
inline

clear function

◆ generateBinaryTree()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::generateBinaryTree ( COMM::ParticleGroup< Tparticle, Tpcm > &  _particles,
const Float  _G 
)
inline

generate binary tree for the particle group

Parameters
[in]_particlesparticle group
[in]_Ggravitational constant

◆ getBinaryTreeRoot()

template<class Tparticle , class Tpcm >
BinaryTree<Tparticle>& AR::Information< Tparticle, Tpcm >::getBinaryTreeRoot ( ) const
inline

get the root of binary tree

Here is the caller graph for this function:

◆ printColumn()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::printColumn ( std::ostream &  _fout,
const int  _width = 20 
)
inline

print data of class members using column style

print data of class members in one line for column style. Notice no newline is printed at the end

Parameters
[out]_foutstd::ostream output object
[in]_widthprint width (defaulted 20)

◆ printColumnTitle()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::printColumnTitle ( std::ostream &  _fout,
const int  _width = 20 
)
inline

print titles of class members using column style

print titles of class members in one line for column style

Parameters
[out]_foutstd::ostream output object
[in]_widthprint width (defaulted 20)

◆ readBinary()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::readBinary ( FILE *  _fin)
inline

read class data to file with binary format

Parameters
[in]_finFILE type file for reading

◆ reserveMem()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::reserveMem ( const int  _nmax)
inline

reserve memory of binarytree list

◆ writeBinary()

template<class Tparticle , class Tpcm >
void AR::Information< Tparticle, Tpcm >::writeBinary ( FILE *  _fout) const
inline

write class data to file with binary format

Parameters
[in]_foutFILE type file for output

Member Data Documentation

◆ binarytree

template<class Tparticle , class Tpcm >
COMM::List<BinaryTree<Tparticle> > AR::Information< Tparticle, Tpcm >::binarytree

fix step option for integration

◆ ds

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::ds

◆ fix_step_option

template<class Tparticle , class Tpcm >
FixStepOption AR::Information< Tparticle, Tpcm >::fix_step_option

◆ r_break_crit

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::r_break_crit

offset of time to obtain real physical time (real time = TimeTransformedSymplecticIntegrator:time_ + info.time_offset)

◆ time_offset

template<class Tparticle , class Tpcm >
Float AR::Information< Tparticle, Tpcm >::time_offset

initial step size for integration


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