HNBody  Version 1.0.10
Data Structures | Macros | Typedefs | Functions | Variables
expert.h File Reference

Low-level function declarations used by the HNBody package. More...

#include <mutils/platform.h>
#include "hnbody/kernel.h"
Include dependency graph for expert.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  hnb_ctrl_struct
 User control callback information. More...
 

Macros

#define vecsize   (3*sizeof(double))
 
#define ACE_MAX   28
 
#define CS_A   1.3512071919596576
 
#define CS_B   -1.7024143839193153
 
#define CS_BA   -1.2599210498948732
 
#define NZP(sys)   (sys->Ltype==LWP_ZWP ? sys->nHL : sys->nH)
 
#define NCM(sys)   (sys->fixedH>0 ? sys->nH : sys->nHL)
 
#define FIXEDH(sys)   (sys->fixedH && sys->nH==sys->nHL)
 
#define FIXEDZP(sys)
 
#define IS_BARY(it)   (((it)&Origin)==Barycentric)
 
#define IS_BODY(it)   (((it)&Origin)==Bodycentric)
 
#define IS_JACOBI(it)   (((it)&Origin)==Jacobi)
 
#define IS_REG(it)   (((it)&Regularized)!=0)
 
#define IS_MOD(it)   (((it)&TWform)!=0)
 
#define IS_4TH(it)   (((it)&Order)!=0)
 
#define IS_PN(it)   (((it)&Gravity)==PostNewtonian)
 
#define IS_CS(it)   (((it)&Order)==Order4)
 
#define IS_KD(it)   (((it)&Splitting)==Kick_Drift)
 
#define IS_DK(it)   (((it)&Splitting)==Drift_Kick)
 
#define IS_SD(it)   (((it)&Splitting)==Shift_Drift)
 
#define IMAP(j, nL, nHL)   ((j)<=(nL) ? ((j) ? (nHL)-(j) : 0) : (j)-(nL))
 
#define ABSLT(x, y)   (fabs(x)<fabs(y))
 
#define ABSLE(x, y)   (fabs(x)<=fabs(y))
 
#define ABSGT(x, y)   (fabs(x)>fabs(y))
 
#define ABSGE(x, y)   (fabs(x)>=fabs(y))
 
#define memzero(x, n)   memset(x, 0, (n)*sizeof(double))
 
#define hnb_ksum(S, C, Xj)   do {double Y=Xj-C, T=S+Y; C=(T-S)-Y; S=T;} while (0)
 

Typedefs

typedef struct hnb_ctrl_struct hnb_ctrl_t
 User control callback information. More...
 

Functions

hnb_vector_t hnb_newvec (hnb_data_t *sys)
 
hnb_vector_t hnb_newvec_unlocked (hnb_data_t *sys)
 
DLLSPEC hnb_data_thnb_raw_init (double tinit, const double xinit[][3], const double vinit[][3], const double m[], const int id[], const int jindex[], int nH, int nL, int nZ, hnb_LWP_t ltype, int fixedH, int enc, int prune, const double renc[], const double rcapt[], double G, double c, double Msun, double J2, double J4, double J6, double obRadius, double dm, double dzH, double dzZ, double eps, hnb_integ_t integ, hnb_integcoord_t itH, hnb_integcoord_t itZ, hnb_extra_t extra_kick, hnb_extra_t extra_shift, hnb_drift_t extra_drift, hnb_derivs_t extra_derivs, hnb_energy_t extra_energy, int nthreads)
 Allocates and initializes a system data structure. More...
 
DLLSPEC void hnb_setup_output (hnb_data_t *sys)
 
DLLSPEC void hnb_drift (hnb_data_t *sys, double dz, int output)
 
DLLSPEC void hnb_kick (hnb_data_t *sys, double dz, int output)
 
DLLSPEC void hnb_shift (hnb_data_t *sys, double dz, int output)
 
DLLSPEC void * hnb_mtsort (void *vsys)
 
DLLSPEC void hnb_init_threads (hnb_data_t *sys, int npart)
 
DLLSPEC void hnb_delvec (hnb_vector_t vec, hnb_data_t *sys)
 
DLLSPEC void hnb_delvec_unlocked (hnb_vector_t vec, hnb_data_t *sys)
 
DLLSPEC void hnb_halfstep (hnb_data_t *sys)
 
DLLSPEC void ode_step (hnb_data_t *sys)
 
DLLSPEC void zwp_iostep (hnb_data_t *, int)
 
DLLSPEC void zwp_step (hnb_data_t *)
 
DLLSPEC void zwp_regdata (double *, double *, double *, const double *, const double *, hnb_data_t *)
 
DLLSPEC void hcoord_to_zcoord (double(*)[3], double(*)[3], hnb_data_t *)
 
DLLSPEC void hl_dkick (double dv[][3], double x[][3], double dt, hnb_data_t *sys)
 
DLLSPEC void hl_drift (hnb_data_t *sys, double dz)
 
DLLSPEC void hl_kick0 (hnb_data_t *sys, double dz)
 
DLLSPEC void hl_shift (hnb_data_t *sys, double dz)
 
DLLSPEC void hl_sync (hnb_data_t *, double)
 
DLLSPEC void hl_reinit (hnb_data_t *, int)
 
DLLSPEC void hl_etas (double *, double *, const double *, const int *, int)
 
DLLSPEC void hl_iostep (hnb_data_t *, int)
 
DLLSPEC void hl_step0 (hnb_data_t *, double, int)
 
DLLSPEC void hl_inter (double[], double[], double, int, hnb_data_t *)
 
DLLSPEC void hnb_normalize (double x[][3], double dx[][3], int n)
 
DLLSPEC void hnb_partition (int *imin, int *imax, int ipart, int npart)
 
DLLSPEC void oblate_kick (double(*dv)[3], double(*x)[3], const double m[], int n, double dt, hnb_data_t *sys)
 
DLLSPEC void postn_bary2pseudo (double(*vpseudo)[3], double(*xbody)[3], double(*vbary)[3], const double m[], int n, double G, double c2)
 
DLLSPEC void postn_pseudo2bary (double(*vbary)[3], double(*xbody)[3], double(*vpseudo)[3], const double m[], int n, double G, double c2)
 
DLLSPEC void postn_jacobi2pseudo (double(*dv)[3], double(*x)[3], double(*v)[3], const double ieta[], const int imap[], int n, double GM, double c2)
 
DLLSPEC void postn_pseudo2jacobi (double(*dv)[3], double(*x)[3], double(*v)[3], const double ieta[], const int imap[], int n, double GM, double c2)
 
DLLSPEC void postn_shift (double(*dx)[3], double(*v)[3], const double m[], int n, double dt, double G, double c2, hnb_integcoord_t it)
 
DLLSPEC void postn_kick (double(*dv)[3], double(*x)[3], const double m[], const double ieta[], const int imap[], int n, double dt, double G, double c2, hnb_integcoord_t it)
 
DLLSPEC double oblate_energy (double(*x)[3], const double m[], int n, hnb_data_t *sys)
 
DLLSPEC double postn_energy (double(*x)[3], double(*v)[3], hnb_data_t *sys)
 
DLLSPEC int hnb_extra_index (const hnb_data_t *)
 
DLLSPEC int create_imap (int imap[], const int jindex[], int nhl)
 
DLLSPEC int hl_filltab (hnb_data_t *, double)
 
DLLSPEC int close_enc (const double xc0[], const double vc0[], const double xc1[], const double vc1[], const double xp0[], const double vp0[], const double xp1[], const double vp1[], double renc, double GM, double dt)
 
DLLSPEC double encke_g (double)
 Computes 1-(1+2*q)^-3/2, maintaining relative precision for |q|<<1.
 

Variables

DLLSPEC hnb_ctrl_t hnb_user_ctrl_data
 Passed to hnb_user_exam_func to allow increased user integration control.
 
DLLSPEC void(* hnb_user_exam_func )(hnb_ctrl_t *ctrl, double t, const double x[][3], const double v[][3], const double m[], int n, const hnb_data_t *sys)
 
DLLSPEC void(*)(*) hnb_user_prune_func (hnb_data_t *sys, int index)
 
DLLSPEC void(*)(*)(*) hnb_user_init_func (hnb_data_t *sys, const char *file)
 
DLLSPEC void(*)(*)(*)(*) hnb_user_exit_func (hnb_data_t *sys, const char *file)
 
DLLSPEC void(*)(*)(*)(*)(*) hnb_user_save_func (hnb_data_t *sys, double key)
 
DLLSPEC void(*)(*)(*)(*)(*)(*) hnb_user_rest_func (hnb_data_t *sys, double key)
 
DLLSPEC void(*)(*)(*)(*)(*)(*)(*) hnb_user_init_step (hnb_data_t *sys)
 
DLLSPEC void(*)(*)(*)(*)(*)(*)(*)(*) hnb_user_half_step (hnb_data_t *sys)
 
DLLSPEC const char * hnb_user_init_file
 

Detailed Description

Id
expert.h,v 1.19 2014/09/03 13:23:38 rauch Exp
Author
Kevin P. Rauch
Warning
These interfaces are for experts only!!!!

Macro Definition Documentation

#define FIXEDZP (   sys)
Value:
( sys->fixedH && \
(sys->nH==sys->nHL || sys->Ltype==LWP_NONE))
Definition: kernel.h:84

Typedef Documentation

typedef struct hnb_ctrl_struct hnb_ctrl_t

This structure is passed to hnb_user_exam_func to allow users to modify actions performed by HNBody during the following step.

The prune vector allows users to control the pruning of individual particles each step. Setting an element to -1, 0, or +1 requests the default action (significantly hyperbolic particles are pruned), disables pruning, or requests immediate pruning, respectively. The contents of this vector is automatically reset to -1 before each call to hnb_user_exam_func. Note that the input file must specify PruneCollisions = true for pruning requests to be honored and that (currently) only Symplectic integrations support pruning.

Warning
prune may be NULL if pruning is not supported.

Function Documentation

DLLSPEC hnb_data_t* hnb_raw_init ( double  tinit,
const double  x[][3],
const double  v[][3],
const double  m[],
const int  id[],
const int  jindex[],
int  nh,
int  nl,
int  nz,
hnb_LWP_t  Ltype,
int  fixedH,
int  enc,
int  prune,
const double  renc[],
const double  rcapt[],
double  G,
double  c,
double  Msun,
double  J2,
double  J4,
double  J6,
double  obRad,
double  dm,
double  dzH,
double  dzZ,
double  eps,
hnb_integ_t  integ,
hnb_integcoord_t  itH,
hnb_integcoord_t  itZ,
hnb_extra_t  extra_kick,
hnb_extra_t  extra_shift,
hnb_drift_t  extra_drift,
hnb_derivs_t  extra_derivs,
hnb_energy_t  extra_energy,
int  nthreads 
)

This routine initializes a dynamically allocated system data structure.

On input, x and v are the bodies' position and velocity in some inertial reference frame at time tinit, and m are their masses. Index 0 (zero) always refers to the central, dominant mass! x[1]..x[nh-1] are the other "big" masses and are best ordered according to distance from x[0]; x[nh]..x[n-1] are "small" masses and can appear in any order, except that truly "massless" bodies should appear last. tH[0]=tZ[0]=0 is set on input and should not be altered. Heap space is allocated for the returned structure and its internal arrays (x, v, m, renc, rcapt, id, and jindex are copied).

Note
All values are assumed to be in a consistent set of units; integrations have no knowledge of the absolute units used.
Parameters
tinitis the initial epoch of integration.
xcontains the initial positions of the bodies.
vcontains the initial velocities of the bodies.
mcontains the masses of the bodies.
idcontains the ID numbers of the bodies.
jindexcontains the Jacobi indices of the bodies, if applicable.
pruneis non-zero iff escaping bodies should be removed from the integration.
nhis the number of HWPs.
nlis the number of LWPs.
nzis the number of ZWPs.
Gis the Newtonian gravitational constant.
cis the speed of light.
extra_kickis an optional extra_kick() term to include in the integration; use NULL if none (Symplectic only).
extra_shiftis an optional extra_shift() term to include in the integration; use NULL if none (Symplectic only).
extra_driftis an optional extra_drift() term to include in the integration; use NULL if none (Symplectic only).
extra_derivsis an optional extra_derivs() term to include in the integration; use NULL if none (ODE only).
extra_energyis an optional extra_energy() term to include in the integration; use NULL if none.
nthreadsis the number of threads across which to partition data calculations (currently non-functional).
Returns
a valid, initialized integration structure, else NULL.

References Barycentric, Bodycentric, Corrected, Corrected2, Corrected4, Drift_Kick, Integrator, Jacobi, Kick_Drift, LWP_NONE, LWP_ZWP, ModBodycentric, Order2, Order4, hnb_ctrl_struct::prune, RegBarycentric, RungeKutta, Shift_Drift, and Symplectic.

Referenced by hnb_opts_init().

Variable Documentation

DLLSPEC void(* hnb_user_exam_func)(hnb_ctrl_t *ctrl, double t, const double x[][3], const double v[][3], const double m[], int n, const hnb_data_t *sys)

Called just prior to each integration step; for the ODE integrators, these steps refer to the smaller, adaptive steps and not the top-level hnb_step().

Referenced by hnb_interface().