![]() |
![]() |
![]() |
PHOEBE Reference Manual | ![]() |
---|---|---|---|---|
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);
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* phoebe_star_surface_new ();
Initializes a new PHOEBE_star_surface structure.
Returns : |
pointer to the newly initialized PHOEBE_star_surface. |
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
).
|
star surface elements to be allocated |
|
latitude raster size (number of latitude circles) |
Returns : |
PHOEBE_error_code. |
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.
|
stellar surface to be rasterized |
|
raster size |
Returns : |
PHOEBE_error_code. |
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 for which to compute the local radii |
|
surface equipotential value |
|
mass ratio |
|
instantaneous separation |
|
asynchronicity parameter |
Returns : |
PHOEBE_error_code. |
int phoebe_star_surface_compute_gradients (PHOEBE_star_surface *surface, double q, double D, double F);
|
surface for which to compute local gradients |
|
mass ratio |
|
instantaneous separation |
|
asynchronicity parameter |
Returns : |
PHOEBE_error_code |
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 for which to compute cos(beta) |
Returns : |
PHOEBE_error_code |
int phoebe_star_surface_compute_cos_coordinates (PHOEBE_star_surface *surface);
|
|
Returns : |
int phoebe_star_surface_compute_pos_coordinates (PHOEBE_star_surface *surface, double incl, double phase);
|
|
|
|
|
|
Returns : |
int phoebe_star_surface_free (PHOEBE_star_surface *surface);
Frees memory allocated for the star surface surface
and its fields.
|
star surface to be freed |
Returns : |
PHOEBE_error_code. |
int phoebe_star_effective_radius (PHOEBE_star *star, double *radius);
|
|
|
|
Returns : |
int phoebe_star_volume (PHOEBE_star *star, double *volume);
|
|
|
|
Returns : |
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.
|
surface equipotential value |
|
instantaneous separation |
|
mass ratio |
Returns : |
polar 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.
|
polar radius |
|
mass ratio |
|
instantaneous separation |
|
asynchronicity parameter |
|
direction cosine, lambda = sin(theta) cos(phi)
|
|
direction cosine, nu = cos(theta)
|
Returns : |
radius in units of semi-major axis. |