![]() |
![]() |
![]() |
PHOEBE Reference Manual | ![]() |
---|---|---|---|---|
enum PHOEBE_cost_function; double frac (double x); int diff (const void *a, const void *b); int diff_int (const void *a, const void *b); int phoebe_interpolate (int N, double *x, double *lo, double *hi, PHOEBE_type type, ...); int phoebe_cf_compute (double *cfval, PHOEBE_cost_function cf, PHOEBE_vector *syndep, PHOEBE_vector *obsdep, PHOEBE_vector *obsweight, double scale); int phoebe_join_chi2 (double *chi2, PHOEBE_vector *chi2s, PHOEBE_vector *weights); #define wd_lc (atmtab,pltab,lcin,request,vertno,L3perc,indeps,deps,poscoy,poscoz,params) int phoebe_compute_lc_using_wd (PHOEBE_curve *curve, PHOEBE_vector *indep, char *lcin); int phoebe_compute_rv1_using_wd (PHOEBE_curve *rv1, PHOEBE_vector *indep, char *lcin); int phoebe_compute_rv2_using_wd (PHOEBE_curve *rv2, PHOEBE_vector *indep, char *lcin); int phoebe_compute_pos_using_wd (PHOEBE_vector *poscoy, PHOEBE_vector *poscoz, char *lcin, double phase); int call_wd_to_get_logg_values (double *logg1, double *logg2); int phoebe_calculate_level_correction (double *alpha, PHOEBE_curve *syn, PHOEBE_curve *obs); int phoebe_calculate_gamma_correction (double *gamma, PHOEBE_curve *syn, PHOEBE_curve *obs); double phoebe_calculate_pot1 (bool ELLIPTIC, double D, double q, double r, double F, double lambda, double nu); double phoebe_calculate_pot2 (bool ELLIPTIC, double D, double q, double r, double F, double lambda, double nu); int phoebe_calculate_masses (double sma, double P, double q, double *M1, double *M2); int phoebe_calculate_critical_potentials (double q, double F, double e, double *L1crit, double *L2crit); int phoebe_compute_critical_phases (double *pp, double *scp, double *icp, double *anp, double *dnp, double perr0, double ecc, double pshift); int calculate_weighted_sum (double *sum, PHOEBE_vector *dep, PHOEBE_vector *weight); int calculate_weighted_average (double *average, PHOEBE_vector *dep, PHOEBE_vector *weight); int calculate_weighted_sigma (double *sigma, PHOEBE_vector *dep, PHOEBE_vector *weight); double intern_calculate_phase_from_ephemeris (double hjd, double hjd0, double period, double dpdt, double pshift); int transform_hjd_to_phase (PHOEBE_vector *vec, double hjd0, double period, double dpdt, double pshift); int transform_phase_to_hjd (PHOEBE_vector *vec, double hjd0, double period, double dpdt, double pshift); int transform_magnitude_to_flux (PHOEBE_vector *vec, double mnorm); int transform_magnitude_sigma_to_flux_sigma (PHOEBE_vector *weights, PHOEBE_vector *fluxes); int transform_flux_to_magnitude (PHOEBE_vector *vec, double mnorm); int transform_flux_sigma_to_magnitude_sigma (PHOEBE_vector *weights, PHOEBE_vector *fluxes); int normalize_kms_to_orbit (PHOEBE_vector *vec, double sma, double period); int transform_sigma_to_weight (PHOEBE_vector *vec); int transform_weight_to_sigma (PHOEBE_vector *vec); int calculate_main_sequence_parameters (double T1, double T2, double P0, double *L1, double *L2, double *M1, double *M2, double *q, double *a, double *R1, double *R2, double *Omega1, double *Omega2); int calculate_synthetic_scatter_seed (double *seed); int apply_extinction_correction (PHOEBE_curve *curve, double A); int apply_third_light_correction (PHOEBE_curve *curve, PHOEBE_el3_units el3units, double el3value); int apply_interstellar_extinction_correction (PHOEBE_vector *wavelength, PHOEBE_vector *spectrum, double R, double E); int calculate_teff_from_bv_index (int star_type, double bv, double *teff); int phoebe_wd_model_from_phoebe_model_parameter (); bool phoebe_phsv_constrained (int wd_model); bool phoebe_pcsv_constrained (int wd_model);
typedef enum PHOEBE_cost_function { PHOEBE_CF_STANDARD_DEVIATION, PHOEBE_CF_WEIGHTED_STANDARD_DEVIATION, PHOEBE_CF_SUM_OF_SQUARES, PHOEBE_CF_EXPECTATION_CHI2, PHOEBE_CF_CHI2 } PHOEBE_cost_function;
int phoebe_interpolate (int N, double *x, double *lo, double *hi, PHOEBE_type type, ...);
This is a general multi-dimensional linear interpolation function. It should be reasonably optimized - I gave it a lot of thought. It should also be reasonably well tested.
The order of nodes and function values is very important; because of optimization, its order is inverse-binary:
lo[0] = lo[par1], lo[1] = lo[par2], ..., lo[N-1] = lo[parN]
hi[0] = hi[par1], hi[1] = hi[par2], ..., hi[N-1] = hi[parN]
fv[0] = fv (0 0 0 ... 0)
fv[1] = fv (1 0 0 ... 0)
fv[2] = fv (0 1 0 ... 0)
fv[3] = fv (1 1 0 ... 0)
fv[4] = fv (0 0 1 ... 0)
fv[5] = fv (1 0 1 ... 0)
fv[6] = fv (0 1 1 ... 0)
fv[7] = fv (1 1 1 ... 0)
.....................
where 0 and 1 are used for lower and upper node values, respectively, listed in a consecutive parameter order.
The function *modifies* the passed array fv[], so make sure you copy its contents if it is needed later. The result of the interpolation is contained in the first element of that array, fv[0].
|
dimension of the interpolation space |
|
N -dimensional vector to the point of interest
|
|
N -dimensional vector of lower node values
|
|
N -dimensional vector of upper node values
|
|
PHOEBE_type od interpolation values |
|
(2^N)-dimensional vector of node function values of type type .
|
Returns : |
PHOEBE_error_code. |
int phoebe_cf_compute (double *cfval, PHOEBE_cost_function cf, PHOEBE_vector *syndep, PHOEBE_vector *obsdep, PHOEBE_vector *obsweight, double scale);
Computes the cost function value cfval
of the passed cost function
cf
. The residuals are computed from vectors syndep
and obsdep
. If
the cost function is weighted, each residual is multiplied by the
inverse square of the individual obsweight
value. Since the residuals
for different curves are usually compared, a scaling constant scale
can be used to renormalize the data. The scale
is usually computed as
4\pi/(L1+L2+4\piL3).
Cost function cf
is of the following:
PHOEBE_CF_STANDARD_DEVIATION
PHOEBE_CF_WEIGHTED_STANDARD_DEVIATION
PHOEBE_CF_SUM_OF_SQUARES
PHOEBE_CF_EXPECTATION_CHI2
PHOEBE_CF_CHI2
|
the computed cost function value. |
|
a PHOEBE_cost_function to be evaluated. |
|
a PHOEBE_vector of the model data. |
|
a PHOEBE_vector of the observed data. |
|
|
|
a scaling constant for computing the residuals. |
Returns : |
a PHOEBE_error_code. |
int phoebe_join_chi2 (double *chi2, PHOEBE_vector *chi2s, PHOEBE_vector *weights);
|
|
|
|
|
|
Returns : |
#define wd_lc(atmtab,pltab,lcin,request,vertno,L3perc,indeps,deps,poscoy,poscoz,params) lc_(atmtab,pltab,lcin,request,vertno,L3perc,indeps,deps,poscoy,poscoz,params,strlen(atmtab),strlen(pltab),strlen(lcin))
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
int phoebe_compute_lc_using_wd (PHOEBE_curve *curve, PHOEBE_vector *indep, char *lcin);
Uses WD's LC code through a FORTRAN wrapper to obtain the fluxes.
PHOEBE_curve curve
must be initialized.
|
PHOEBE_curve placeholder for computed fluxes |
|
an array of independent variable values (HJDs or phases) |
|
WD's lci filename |
Returns : |
PHOEBE_error_code. |
int phoebe_compute_rv1_using_wd (PHOEBE_curve *rv1, PHOEBE_vector *indep, char *lcin);
Uses WD's LC code through a FORTRAN wrapper to obtain the primary
star radial velocities. PHOEBE_curve rv1
must be initialized.
|
PHOEBE_curve placeholder for computed radial velocities |
|
an array of independent variable values (HJDs or phases); |
|
|
Returns : |
PHOEBE_error_code. |
int phoebe_compute_rv2_using_wd (PHOEBE_curve *rv2, PHOEBE_vector *indep, char *lcin);
|
|
|
|
|
|
Returns : |
int phoebe_compute_pos_using_wd (PHOEBE_vector *poscoy, PHOEBE_vector *poscoz, char *lcin, double phase);
Uses WD's LC code through a FORTRAN wrapper to obtain the plane-of-sky
coordinates. PHOEBE_vectors poscoy
and poscoz
must be initialized.
|
PHOEBE_vector placeholder for the plane-of-sky y coordinates |
|
PHOEBE_vector placeholder for the plane-of-sky z coordinates |
|
|
|
phase at which the plane-of-sky coordinates are computed |
Returns : |
PHOEBE_error_code. |
int call_wd_to_get_logg_values (double *logg1, double *logg2);
|
|
|
|
Returns : |
int phoebe_calculate_level_correction (double *alpha, PHOEBE_curve *syn, PHOEBE_curve *obs);
Computes the correction alpha
by solving:
\sum_i w_i (o_i - alpha
c_i)^2 = min
where o_i are observed data points with individual weights w_i, and c_i are the computed fluxes with the original value of L1. Solving this system requires solving the following least squares equation:
alpha
= (A W A^T)^{-1} A W b = \sum_i w_i o_i c_i / \sum_i w_i o_i^2
|
placeholder for level correction: L1 -> L1/alpha , L2 -> L2/alpha
|
|
synthetic (model) light curve |
|
observed light curve |
Returns : |
PHOEBE_error_code. |
int phoebe_calculate_gamma_correction (double *gamma, PHOEBE_curve *syn, PHOEBE_curve *obs);
Computes the correction gamma
by solving:
\sum_i w_i (o_i - c_i - gamma
)^2 = min
where o_i are observed data points with individual weights w_i, and c_i
are the computed fluxes with the original value of gamma
. Solving this
system requires solving the following least squares equation:
gamma
= (A W A^T)^{-1} A W b = \sum_i w_i o_i / \sum_i w_i
|
placeholder for gamma (center-of-mass) velocity correction: v -> v + gamma |
|
synthetic (model) light curve |
|
observed light curve |
Returns : |
PHOEBE_error_code. |
double phoebe_calculate_pot1 (bool ELLIPTIC, double D, double q, double r, double F, double lambda, double nu);
Calculates the primary star surface potential according to Eq.(3.16) in PHOEBE scientific reference.
|
elliptic orbit switch: 1 for elliptic orbits, 0 for circular |
|
instantaneous separation in units of semi-major axis |
|
mass ratio |
|
effective radius of the star |
|
asynchronicity parameter |
|
direction cosine, lambda = sin(theta) cos(phi)
|
|
direction cosine, nu = cos(theta)
|
Returns : |
value of the primary star surface potential |
double phoebe_calculate_pot2 (bool ELLIPTIC, double D, double q, double r, double F, double lambda, double nu);
Calculates the secondary star surface potential by transforming the coordinate system to the center of the secondary star according to Eq.(3.15) in PHOEBE scientific reference, and using Eq.(3.16) to do the calculation.
|
elliptic orbit switch: 1 for elliptic orbits, 0 for circular |
|
instantaneous separation in units of semi-major axis |
|
mass ratio |
|
effective radius of the star |
|
asynchronicity parameter |
|
direction cosine, lambda = sin(theta) cos(phi)
|
|
direction cosine, nu = cos(theta)
|
Returns : |
value of the secondary star surface potential |
int phoebe_calculate_masses (double sma, double P, double q, double *M1, double *M2);
Computes the mass of the primary and the secondary star, in units of solar mass.
|
semi-major axis in solar radii |
|
orbital period in days |
|
mass ratio, M2/M1 |
|
placeholder for the primary star mass |
|
placeholder for the secondary star mass |
Returns : |
PHOEBE_error_code. |
int phoebe_calculate_critical_potentials (double q, double F, double e, double *L1crit, double *L2crit);
Calculates the value of the primary surface potential in Lagrange points L1 and L2.
|
mass ratio |
|
asynchronicity parameter |
|
eccentricity |
|
placeholder for the inner (Roche) lobe value of the potential |
|
placeholder for the outer lobe value of the potential |
Returns : |
PHOEBE_error_code |
int phoebe_compute_critical_phases (double *pp, double *scp, double *icp, double *anp, double *dnp, double perr0, double ecc, double pshift);
Computes critical phases by solving the Kepler problem for critical true anomalies:
Periastron passage: pp
= perr0
/(2pi) - 0.25 + pshift
Superior conjunction: T = pi/2-perr0
Inferior conjunction: T = 3pi/2-perr0
Ascending node: T = -perr0
Descending node: T = pi-perr0
|
placeholder for periastron passage phase |
|
placeholder for superior conjunction phase |
|
placeholder for inferior conjunction phase |
|
placeholder for ascending node phase |
|
placeholder for descending node phase |
|
argument of periastron, in radians |
|
orbital eccentricity |
|
|
Returns : |
PHOEBE_error_code. |
int calculate_weighted_sum (double *sum, PHOEBE_vector *dep, PHOEBE_vector *weight);
|
|
|
|
|
|
Returns : |
int calculate_weighted_average (double *average, PHOEBE_vector *dep, PHOEBE_vector *weight);
|
|
|
|
|
|
Returns : |
int calculate_weighted_sigma (double *sigma, PHOEBE_vector *dep, PHOEBE_vector *weight);
|
|
|
|
|
|
Returns : |
double intern_calculate_phase_from_ephemeris (double hjd, double hjd0, double period, double dpdt, double pshift);
|
|
|
|
|
|
|
|
|
|
Returns : |
int transform_hjd_to_phase (PHOEBE_vector *vec, double hjd0, double period, double dpdt, double pshift);
|
|
|
|
|
|
|
|
|
|
Returns : |
int transform_phase_to_hjd (PHOEBE_vector *vec, double hjd0, double period, double dpdt, double pshift);
|
|
|
|
|
|
|
|
|
|
Returns : |
int transform_magnitude_to_flux (PHOEBE_vector *vec, double mnorm);
|
|
|
|
Returns : |
int transform_magnitude_sigma_to_flux_sigma (PHOEBE_vector *weights, PHOEBE_vector *fluxes);
|
|
|
|
Returns : |
int transform_flux_to_magnitude (PHOEBE_vector *vec, double mnorm);
|
|
|
|
Returns : |
int transform_flux_sigma_to_magnitude_sigma (PHOEBE_vector *weights, PHOEBE_vector *fluxes);
|
|
|
|
Returns : |
int normalize_kms_to_orbit (PHOEBE_vector *vec, double sma, double period);
|
|
|
|
|
|
Returns : |
int calculate_main_sequence_parameters (double T1, double T2, double P0, double *L1, double *L2, double *M1, double *M2, double *q, double *a, double *R1, double *R2, double *Omega1, double *Omega2);
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Returns : |
int calculate_synthetic_scatter_seed (double *seed);
|
|
Returns : |
int apply_extinction_correction (PHOEBE_curve *curve, double A);
|
|
|
|
Returns : |
int apply_third_light_correction (PHOEBE_curve *curve, PHOEBE_el3_units el3units, double el3value);
|
|
|
|
|
|
Returns : |
int apply_interstellar_extinction_correction (PHOEBE_vector *wavelength, PHOEBE_vector *spectrum, double R, double E);
|
|
|
|
|
|
|
|
Returns : |
int calculate_teff_from_bv_index (int star_type, double bv, double *teff);
|
|
|
|
|
|
Returns : |
int phoebe_wd_model_from_phoebe_model_parameter ();
Translates the parameter "phoebe_model" into the WD model number (-1 to 6).
Returns : |
wd_model .
|
bool phoebe_phsv_constrained (int wd_model);
Indicates whether the primary potential is constrained by the choice of model.
|
Model number used by WD (-1 to 6) |
Returns : |
TRUE/FALSE. |
bool phoebe_pcsv_constrained (int wd_model);
Indicates whether the secondary potential is constrained by the choice of model.
|
Model number used by WD (-1 to 6) |
Returns : |
TRUE/FALSE. |