ASCOT5
|
Library of Ascot5 functions for external use. More...
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <hdf5.h>
#include <math.h>
#include "ascot5.h"
#include "gitver.h"
#include "simulate.h"
#include "B_field.h"
#include "E_field.h"
#include "plasma.h"
#include "wall.h"
#include "neutral.h"
#include "boozer.h"
#include "mhd.h"
#include "asigma.h"
#include "consts.h"
#include "physlib.h"
#include "simulate/mccc/mccc_coefs.h"
#include "hdf5_interface.h"
#include "hdf5io/hdf5_helpers.h"
#include "hdf5io/hdf5_bfield.h"
#include "hdf5io/hdf5_efield.h"
#include "hdf5io/hdf5_plasma.h"
#include "hdf5io/hdf5_wall.h"
#include "hdf5io/hdf5_neutral.h"
#include "hdf5io/hdf5_boozer.h"
#include "hdf5io/hdf5_mhd.h"
Go to the source code of this file.
Functions | |
void | libascot_B_field_eval_B_dB (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *BR, real *Bphi, real *Bz, real *BR_dR, real *BR_dphi, real *BR_dz, real *Bphi_dR, real *Bphi_dphi, real *Bphi_dz, real *Bz_dR, real *Bz_dphi, real *Bz_dz) |
Evaluate magnetic field vector and derivatives at given coordinates. | |
void | libascot_B_field_eval_rho (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *rho, real *drhodpsi, real *psi, real *dpsidr, real *dpsidphi, real *dpsidz) |
Evaluate normalized poloidal flux at given coordinates. | |
void | libascot_B_field_get_axis (sim_data *sim, int Neval, real *phi, real *Raxis, real *zaxis) |
Get magnetic axis at given coordinates. | |
void | libascot_B_field_rhotheta2rz (sim_data *sim, int Neval, real *rho, real *theta, real *phi, real t, int maxiter, real tol, real *r, real *z) |
Map (rho, theta, phi) to (R,z) coordinates. | |
void | libascot_B_field_gradient_descent (sim_data *sim, real psi[1], real rz[2], real step, real tol, int maxiter, int ascent) |
Find psi on axis using the gradient descent method. | |
void | libascot_B_field_gradient_descent_3d (sim_data *sim, real psi[1], real rzphi[3], real phimin, real phimax, real step, real tol, int maxiter, int ascent) |
Find one psi minimum using the gradient descent method for 3D field inside a sector (phimin < phi < phimax). Not guaranteed to find the global minimum inside the given sector. | |
void | libascot_E_field_eval_E (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *ER, real *Ephi, real *Ez) |
Evaluate electric field vector at given coordinates. | |
int | libascot_plasma_get_n_species (sim_data *sim) |
Get number of plasma species. | |
void | libascot_plasma_get_species_mass_and_charge (sim_data *sim, real *mass, real *charge, int *anum, int *znum) |
Get mass and charge of all plasma species. | |
void | libascot_plasma_eval_background (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *dens, real *temp) |
Evaluate plasma density and temperature at given coordinates. | |
void | libascot_neutral_eval_density (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *dens) |
Evaluate neutral density at given coordinates. | |
void | libascot_boozer_eval_psithetazeta (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *psi, real *theta, real *zeta, real *dpsidr, real *dpsidphi, real *dpsidz, real *dthetadr, real *dthetadphi, real *dthetadz, real *dzetadr, real *dzetadphi, real *dzetadz, real *rho) |
Evaluate boozer coordinates and derivatives. | |
void | libascot_boozer_eval_fun (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, real *qprof, real *jac, real *jacB2) |
Evaluate boozer coordinates related quantities. | |
int | libascot_mhd_get_n_modes (sim_data *sim) |
Get number of MHD modes. | |
void | libascot_mhd_get_mode_specs (sim_data *sim, int *nmode, int *mmode, real *amplitude, real *omega, real *phase) |
Get MHD mode amplitude, frequency, phase, and mode numbers. | |
void | libascot_mhd_eval (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, int includemode, real *alpha, real *dadr, real *dadphi, real *dadz, real *dadt, real *Phi, real *dPhidr, real *dPhidphi, real *dPhidz, real *dPhidt) |
Evaluate MHD perturbation potentials. | |
void | libascot_mhd_eval_perturbation (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, int includemode, real *mhd_br, real *mhd_bphi, real *mhd_bz, real *mhd_er, real *mhd_ephi, real *mhd_ez, real *mhd_phi) |
Evaluate MHD perturbation EM-field components. | |
void | libascot_eval_collcoefs (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, int Nv, real *va, real ma, real qa, real *F, real *Dpara, real *Dperp, real *K, real *nu, real *Q, real *dQ, real *dDpara, real *clog, real *mu0, real *mu1, real *dmu0) |
Evaluate collision coefficients. | |
void | libascot_eval_ratecoeff (sim_data *sim, int Neval, real *R, real *phi, real *z, real *t, int Nv, real *va, int Aa, int Za, real ma, int reac_type, real *ratecoeff) |
Evaluate atomic reaction rate coefficient. | |
Library of Ascot5 functions for external use.
Functions in this file allows to evaluate input data and quantities using the same methods as is used in actual simulation.
Definition in file libascot.c.
void libascot_B_field_eval_B_dB | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | BR, | ||
real * | Bphi, | ||
real * | Bz, | ||
real * | BR_dR, | ||
real * | BR_dphi, | ||
real * | BR_dz, | ||
real * | Bphi_dR, | ||
real * | Bphi_dphi, | ||
real * | Bphi_dz, | ||
real * | Bz_dR, | ||
real * | Bz_dphi, | ||
real * | Bz_dz ) |
Evaluate magnetic field vector and derivatives at given coordinates.
sim | simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
BR | output array [T]. |
Bphi | output array [T]. |
Bz | output array [T]. |
BR_dR | output array [T]. |
BR_dphi | output array [T]. |
BR_dz | output array [T]. |
Bphi_dR | output array [T]. |
Bphi_dphi | output array [T]. |
Bphi_dz | output array [T]. |
Bz_dR | output array [T]. |
Bz_dphi | output array [T]. |
Bz_dz | output array [T]. |
Definition at line 65 of file libascot.c.
void libascot_B_field_eval_rho | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | rho, | ||
real * | drhodpsi, | ||
real * | psi, | ||
real * | dpsidr, | ||
real * | dpsidphi, | ||
real * | dpsidz ) |
Evaluate normalized poloidal flux at given coordinates.
sim | simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
rho | output array for the normalized poloidal flux. |
drhodpsi | output array for normalized poloidal flux psi derivative. |
psi | output array for the poloidal flux [Wb]. |
dpsidr | output array for the poloidal flux R derivative [Wb/m]. |
dpsidphi | output array for the poloidal flux phi derivative [Wb/rad]. |
dpsidz | output array for the poloidal flux z derivative [Wb/m]. |
Definition at line 108 of file libascot.c.
void libascot_B_field_get_axis | ( | sim_data * | sim, |
int | Neval, | ||
real * | phi, | ||
real * | Raxis, | ||
real * | zaxis ) |
Get magnetic axis at given coordinates.
sim | simulation data struct |
Neval | number of evaluation points. |
phi | phi coordinates of the evaluation points [rad]. |
Raxis | output array for axis R coordinates. |
zaxis | output array for axis z coordinates. |
Definition at line 141 of file libascot.c.
void libascot_B_field_rhotheta2rz | ( | sim_data * | sim, |
int | Neval, | ||
real * | rho, | ||
real * | theta, | ||
real * | phi, | ||
real | t, | ||
int | maxiter, | ||
real | tol, | ||
real * | r, | ||
real * | z ) |
Map (rho, theta, phi) to (R,z) coordinates.
This function implements the Newton method. If the function fails on a given position, the corresponding (R,z) values in the output arrays are not altered.
sim_data | initialized simulation data struct |
Neval | number of query points. |
rho | the square root of the normalized poloidal flux values. |
theta | poloidal angles [rad]. |
phi | toroidal angles [rad]. |
t | time coordinate (same for all) [s]. |
maxiter | maximum number of iterations in Newton algorithm. |
tol | algorithm is stopped when |rho - rho(r,z)| < tol |
r | output array for R coordinates [m]. |
z | output array for z coordinates [m]. |
Definition at line 173 of file libascot.c.
void libascot_B_field_gradient_descent | ( | sim_data * | sim, |
real | psi[1], | ||
real | rz[2], | ||
real | step, | ||
real | tol, | ||
int | maxiter, | ||
int | ascent ) |
Find psi on axis using the gradient descent method.
Note that the psi value is not returned in case this algorithm fails.
sim_data | initialized simulation data struct |
psi | value of psi on axis if this function did not fail |
rz | initial (R,z) position where also the result is stored |
step | the step size |
tol | the current position is accepted if the distance (in meters) between this and the previous point is below this value |
maxiter | maximum number of iterations before failure |
ascent | if true the algorithm instead ascends to find psi0 (> psi1) |
Definition at line 237 of file libascot.c.
void libascot_B_field_gradient_descent_3d | ( | sim_data * | sim, |
real | psi[1], | ||
real | rzphi[3], | ||
real | phimin, | ||
real | phimax, | ||
real | step, | ||
real | tol, | ||
int | maxiter, | ||
int | ascent ) |
Find one psi minimum using the gradient descent method for 3D field inside a sector (phimin < phi < phimax). Not guaranteed to find the global minimum inside the given sector.
sim_offload_data | initialized simulation offload data struct |
B_offload_array | initialized magnetic field offload data |
psi | value of psi on axis if this function did not fail |
rzphi | initial (R,z,phi) position where also the result is stored |
step | the step size |
tol | the current position is accepted if the distance (in meters) between this and the previous point is below this value |
maxiter | maximum number of iterations before failure |
ascent | if true the algorithm instead ascends to find psi0 (> psi1) |
Definition at line 298 of file libascot.c.
void libascot_E_field_eval_E | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | ER, | ||
real * | Ephi, | ||
real * | Ez ) |
Evaluate electric field vector at given coordinates.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
ER | output array [V/m]. |
Ephi | output array [V/m]. |
Ez | output array [V/m]. |
Definition at line 375 of file libascot.c.
int libascot_plasma_get_n_species | ( | sim_data * | sim | ) |
Get number of plasma species.
param sim_data initialized simulation data struct
Definition at line 399 of file libascot.c.
void libascot_plasma_get_species_mass_and_charge | ( | sim_data * | sim, |
real * | mass, | ||
real * | charge, | ||
int * | anum, | ||
int * | znum ) |
Get mass and charge of all plasma species.
sim_data | initialized simulation data struct |
mass | mass output array [kg]. |
charge | charge output array [C]. |
anum | atomic mass number output array [1]. |
znum | charge number output array [1]. |
Definition at line 412 of file libascot.c.
void libascot_plasma_eval_background | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | dens, | ||
real * | temp ) |
Evaluate plasma density and temperature at given coordinates.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
dens | output array [m^-3]. |
temp | output array [eV]. |
Definition at line 444 of file libascot.c.
void libascot_neutral_eval_density | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | dens ) |
Evaluate neutral density at given coordinates.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
dens | output array [m^-3]. |
Definition at line 481 of file libascot.c.
void libascot_boozer_eval_psithetazeta | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | psi, | ||
real * | theta, | ||
real * | zeta, | ||
real * | dpsidr, | ||
real * | dpsidphi, | ||
real * | dpsidz, | ||
real * | dthetadr, | ||
real * | dthetadphi, | ||
real * | dthetadz, | ||
real * | dzetadr, | ||
real * | dzetadphi, | ||
real * | dzetadz, | ||
real * | rho ) |
Evaluate boozer coordinates and derivatives.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
psi | output array |
theta | output array |
zeta | output array |
dpsidr | output array |
dpsidphi | output array |
dpsidz | output array |
dthetadr | output array |
dthetadphi | output array |
dthetadz | output array |
dzetadr | output array |
dzetadphi | output array |
dzetadz | output array |
rho | output array |
Definition at line 524 of file libascot.c.
void libascot_boozer_eval_fun | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
real * | qprof, | ||
real * | jac, | ||
real * | jacB2 ) |
Evaluate boozer coordinates related quantities.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
qprof | array for storing the (flux averaged) safety factor. |
jac | array for storing the coordinate Jacobian. |
jacB2 | array for storing the coordinate Jacobian multiplied with B^2. |
Definition at line 574 of file libascot.c.
int libascot_mhd_get_n_modes | ( | sim_data * | sim | ) |
Get number of MHD modes.
param sim_data initialized simulation data struct
Definition at line 623 of file libascot.c.
void libascot_mhd_get_mode_specs | ( | sim_data * | sim, |
int * | nmode, | ||
int * | mmode, | ||
real * | amplitude, | ||
real * | omega, | ||
real * | phase ) |
Get MHD mode amplitude, frequency, phase, and mode numbers.
sim_data | initialized simulation data struct |
nmode | output array for toroidal mode number |
mmode | output array for poloidal mode number |
amplitude | output array for mode amplitude |
omega | output array for mode frequency |
phase | output array for mode phase |
Definition at line 638 of file libascot.c.
void libascot_mhd_eval | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
int | includemode, | ||
real * | alpha, | ||
real * | dadr, | ||
real * | dadphi, | ||
real * | dadz, | ||
real * | dadt, | ||
real * | Phi, | ||
real * | dPhidr, | ||
real * | dPhidphi, | ||
real * | dPhidz, | ||
real * | dPhidt ) |
Evaluate MHD perturbation potentials.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
includemode | mode index to include or MHD_INCLUDE_ALL |
alpha | output array |
dadr | output array |
dadphi | output array |
dadz | output array |
dadt | output array |
Phi | output array |
dPhidr | output array |
dPhidphi | output array |
dPhidz | output array |
dPhidt | output array |
Definition at line 678 of file libascot.c.
void libascot_mhd_eval_perturbation | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
int | includemode, | ||
real * | mhd_br, | ||
real * | mhd_bphi, | ||
real * | mhd_bz, | ||
real * | mhd_er, | ||
real * | mhd_ephi, | ||
real * | mhd_ez, | ||
real * | mhd_phi ) |
Evaluate MHD perturbation EM-field components.
sim_data | initialized simulation data struct |
Neval | number of evaluation points. |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
includemode | mode index to include or MHD_INCLUDE_ALL |
mhd_br | output array |
mhd_bphi | output array |
mhd_bz | output array |
mhd_er | output array |
mhd_ephi | output array |
mhd_ez | output array |
mhd_phi | output array |
Definition at line 722 of file libascot.c.
void libascot_eval_collcoefs | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
int | Nv, | ||
real * | va, | ||
real | ma, | ||
real | qa, | ||
real * | F, | ||
real * | Dpara, | ||
real * | Dperp, | ||
real * | K, | ||
real * | nu, | ||
real * | Q, | ||
real * | dQ, | ||
real * | dDpara, | ||
real * | clog, | ||
real * | mu0, | ||
real * | mu1, | ||
real * | dmu0 ) |
Evaluate collision coefficients.
sim_data | initialized simulation data struct |
Neval | number of evaluation points |
R | R coordinates of the evaluation points [m] |
phi | phi coordinates of the evaluation points [rad] |
z | z coordinates of the evaluation points [m] |
t | time coordinates of the evaluation points [s] |
Nv | number of velocity grid points |
va | test particle velocities at the evaluation point [m/s] |
ma | test particle mass [kg] |
qa | test particle charge [C] |
F | output array |
Dpara | output array |
Dperp | output array |
K | output array |
nu | output array |
Q | output array |
dQ | output array |
dDpara | output array |
clog | output array |
mu0 | output array |
mu1 | output array |
dmu0 | output array |
Definition at line 773 of file libascot.c.
void libascot_eval_ratecoeff | ( | sim_data * | sim, |
int | Neval, | ||
real * | R, | ||
real * | phi, | ||
real * | z, | ||
real * | t, | ||
int | Nv, | ||
real * | va, | ||
int | Aa, | ||
int | Za, | ||
real | ma, | ||
int | reac_type, | ||
real * | ratecoeff ) |
Evaluate atomic reaction rate coefficient.
sim_data | initialized simulation data struct |
Neval | number of evaluation points in (R, phi, z, t). |
R | R coordinates of the evaluation points [m]. |
phi | phi coordinates of the evaluation points [rad]. |
z | z coordinates of the evaluation points [m]. |
t | time coordinates of the evaluation points [s]. |
Nv | number of evaluation points in velocity space. |
va | test particle velocities [m/s]. |
Aa | test particle atomic mass number. |
Za | test particle charge number. |
ma | test particle mass. |
reac_type | reaction type |
ratecoeff | output array where evaluated values are stored [1/m^2]. |
Definition at line 872 of file libascot.c.