phoebe_model

phoebe_model

Synopsis

typedef             PHOEBE_star_id;
                    PHOEBE_star_surface;
PHOEBE_star_surface* phoebe_star_surface_new            ();
int                 phoebe_star_surface_alloc           (PHOEBE_star_surface *surface,
                                                         int lat_raster);
int                 phoebe_star_surface_rasterize       (PHOEBE_star_surface *surface,
                                                         int lat_raster);
int                 phoebe_star_surface_compute_radii   (PHOEBE_star_surface *surface,
                                                         double Omega,
                                                         double q,
                                                         double D,
                                                         double F);
int                 phoebe_star_surface_compute_gradients
                                                        (PHOEBE_star_surface *surface,
                                                         double q,
                                                         double D,
                                                         double F);
int                 phoebe_star_surface_compute_cosbeta (PHOEBE_star_surface *surface);
int                 phoebe_star_surface_compute_cos_coordinates
                                                        (PHOEBE_star_surface *surface);
int                 phoebe_star_surface_compute_pos_coordinates
                                                        (PHOEBE_star_surface *surface,
                                                         double incl,
                                                         double phase);
int                 phoebe_star_surface_free            (PHOEBE_star_surface *surface);
                    PHOEBE_star;
int                 phoebe_star_effective_radius        (PHOEBE_star *star,
                                                         double *radius);
int                 phoebe_star_area                    (PHOEBE_star *star,
                                                         double *area);
int                 phoebe_star_volume                  (PHOEBE_star *star,
                                                         double *volume);
double              phoebe_compute_polar_radius         (double Omega,
                                                         double D,
                                                         double q);
double              phoebe_compute_radius               (double rp,
                                                         double q,
                                                         double D,
                                                         double F,
                                                         double lambda,
                                                         double nu);

Description

Details

PHOEBE_star_id

typedef int PHOEBE_star_id;


PHOEBE_star_surface

typedef struct {
	int     elemno;
	double *theta;
	double *phi;
	double *dtheta;
	double *dphi;
	double *rho;
	double *cosbeta;

	struct {
		double x;
		double y;
		double z;
	} *grad;

	/* Center-of-star coordinates: */
	struct {
		double x;
		double y;
		double z;
	} *cos;

	/* Plane-of-sky coordinates: */
	struct {
		double u;
		double v;
		double w;
	} *pos;

	/* These are WD-compatible arrays and may be removed in future: */
	int    *mmsave;
	double *sinth;
	double *costh;
	double *sinphi;
	double *cosphi;
} PHOEBE_star_surface;


phoebe_star_surface_new ()

PHOEBE_star_surface* phoebe_star_surface_new            ();

Initializes a new PHOEBE_star_surface structure.

Returns :

pointer to the newly initialized PHOEBE_star_surface.

phoebe_star_surface_alloc ()

int                 phoebe_star_surface_alloc           (PHOEBE_star_surface *surface,
                                                         int lat_raster);

Computes the number of elements from the passed latitude raster size lat_raster, and allocates memory for local surface quantities. The number of elements on one half of the stellar hemisphere is computed by dividing each of theta = pi/(2lat_raster) latitude circles into 1.3 lat_raster sin(theta) longitude steps and adding them up. Factor 1.3 is introduced to make rasterization along longitude more dense, needed for a better edge determination. This number is then multiplied by 4 to get the whole star surface.

Allocated local surface quantities are: stellar co-latitude (theta), longitude (phi), local radius (rho), and local gradient (grad).

surface :

star surface elements to be allocated

lat_raster :

latitude raster size (number of latitude circles)

Returns :

PHOEBE_error_code.

phoebe_star_surface_rasterize ()

int                 phoebe_star_surface_rasterize       (PHOEBE_star_surface *surface,
                                                         int lat_raster);

Rasterizes surface onto a gridded mesh with lat_raster elements in co- latitude and 1.3 lat_raster sin(theta) elements in longitude. Memory for the raster is allocated by the function and should be freed by the user after use.

surface :

stellar surface to be rasterized

lat_raster :

raster size

Returns :

PHOEBE_error_code.

phoebe_star_surface_compute_radii ()

int                 phoebe_star_surface_compute_radii   (PHOEBE_star_surface *surface,
                                                         double Omega,
                                                         double q,
                                                         double D,
                                                         double F);

Computes local radii of all surface elements of the PHOEBE_star_surface surface. Memory for surface needs to be allocated prior to calling this function.

surface :

surface for which to compute the local radii

Omega :

surface equipotential value

q :

mass ratio

D :

instantaneous separation

F :

asynchronicity parameter

Returns :

PHOEBE_error_code.

phoebe_star_surface_compute_gradients ()

int                 phoebe_star_surface_compute_gradients
                                                        (PHOEBE_star_surface *surface,
                                                         double q,
                                                         double D,
                                                         double F);

surface :

surface for which to compute local gradients

q :

mass ratio

D :

instantaneous separation

F :

asynchronicity parameter

Returns :

PHOEBE_error_code

phoebe_star_surface_compute_cosbeta ()

int                 phoebe_star_surface_compute_cosbeta (PHOEBE_star_surface *surface);

Computes the cosine of the angle between surface normal (normalized gradient of the surface potential) and radius vector. This angle is directly related to surface deformation and thus surface element. The function requires gradients to be computed already.

surface :

surface for which to compute cos(beta)

Returns :

PHOEBE_error_code

phoebe_star_surface_compute_cos_coordinates ()

int                 phoebe_star_surface_compute_cos_coordinates
                                                        (PHOEBE_star_surface *surface);

surface :

Returns :


phoebe_star_surface_compute_pos_coordinates ()

int                 phoebe_star_surface_compute_pos_coordinates
                                                        (PHOEBE_star_surface *surface,
                                                         double incl,
                                                         double phase);

surface :

incl :

phase :

Returns :


phoebe_star_surface_free ()

int                 phoebe_star_surface_free            (PHOEBE_star_surface *surface);

Frees memory allocated for the star surface surface and its fields.

surface :

star surface to be freed

Returns :

PHOEBE_error_code.

PHOEBE_star

typedef struct {
	PHOEBE_star_id       id;
	PHOEBE_star_surface *surface;
} PHOEBE_star;


phoebe_star_effective_radius ()

int                 phoebe_star_effective_radius        (PHOEBE_star *star,
                                                         double *radius);

star :

radius :

Returns :


phoebe_star_area ()

int                 phoebe_star_area                    (PHOEBE_star *star,
                                                         double *area);

star :

area :

Returns :


phoebe_star_volume ()

int                 phoebe_star_volume                  (PHOEBE_star *star,
                                                         double *volume);

star :

volume :

Returns :


phoebe_compute_polar_radius ()

double              phoebe_compute_polar_radius         (double Omega,
                                                         double D,
                                                         double q);

Given the value of surface potential Omega, instantaneous separation D, and mass ratio q, this function computes the value of polar radius in units of semi-major axis. The computation is done iteratively, using the Newton-Raphson scheme.

Omega :

surface equipotential value

D :

instantaneous separation

q :

mass ratio

Returns :

polar radius.

phoebe_compute_radius ()

double              phoebe_compute_radius               (double rp,
                                                         double q,
                                                         double D,
                                                         double F,
                                                         double lambda,
                                                         double nu);

Computes the fractional radius at the given direction iteratively, using the Newton-Raphson scheme.

rp :

polar radius

q :

mass ratio

D :

instantaneous separation

F :

asynchronicity parameter

lambda :

direction cosine, lambda = sin(theta) cos(phi)

nu :

direction cosine, nu = cos(theta)

Returns :

radius in units of semi-major axis.