PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Numerical utilities related to the solution vectors. More...
Interfaces and Types | |
interface | calc_xuq |
Calculates the normal ( \(\cdot_n\)) or geodesic ( \(\cdot_g\)) components of the plasma perturbation \(\vec{\xi}\) or the magnetic perturbation \(\vec{Q} = \nabla \times \left(\vec{x} \times \vec{B}\right)\). More... | |
Functions/Subroutines | |
integer function, public | calc_tot_sol_vec (mds, grid_sol, sol_vec_loc, sol_vec_tot, deriv) |
Calculate the total version of the solution vector from the local version. More... | |
integer function, public | calc_loc_sol_vec (mds, i_min, sol_vec_tot, sol_vec_loc) |
Calculate the local version of the solution vector from the total version. More... | |
integer function, public | interp_v_spline (V_i, V_o, r_i, r_o, extrap, ivs_stat) |
Interpolation for a quantity V using splines. More... | |
Variables | |
logical, public | debug_calc_xuq_arr = .false. |
plot debug information for calc_XUQ_arr More... | |
logical, public | debug_interp_v_spline = .false. |
debug information for interp_v_spline More... | |
Numerical utilities related to the solution vectors.
integer function, public sol_utilities::calc_loc_sol_vec | ( | type(modes_type), intent(in) | mds, |
integer, intent(in) | i_min, | ||
complex(dp), dimension(:,:), intent(in) | sol_vec_tot, | ||
complex(dp), dimension(:,:), intent(inout), allocatable | sol_vec_loc | ||
) |
Calculate the local version of the solution vector from the total version.
[in] | mds | general modes variables |
[in] | i_min | i_min of grid in which variables are tabulated |
[in] | sol_vec_tot | solution vector for all possible resonating modes |
[in,out] | sol_vec_loc | solution vector for local resonating modes |
Definition at line 722 of file sol_utilities.f90.
integer function, public sol_utilities::calc_tot_sol_vec | ( | type(modes_type), intent(in) | mds, |
type(grid_type), intent(in) | grid_sol, | ||
complex(dp), dimension(:,:), intent(in) | sol_vec_loc, | ||
complex(dp), dimension(:,:), intent(inout), allocatable | sol_vec_tot, | ||
integer, intent(in), optional | deriv | ||
) |
Calculate the total version of the solution vector from the local version.
This is the solution vector for all of the possible mode numbers, which can be different from the local mode numbers for X style 2 (fast):
n_mod_X
min_sec_X
and max_sec_X
.X_style
2, as the mode numbers depend on the normal coordinate.If the output variable is not allocated, it is done here.
Optionally a derivative can be requested, depending on the normal discretization order.
[in] | mds | general modes variables |
[in] | grid_sol | solution grid |
[in] | sol_vec_loc | solution vector for local resonating modes |
[in,out] | sol_vec_tot | solution vector for all possible resonating modes |
[in] | deriv | order of derivative |
Definition at line 581 of file sol_utilities.f90.
integer function, public sol_utilities::interp_v_spline | ( | complex(dp), dimension(:,:,:), intent(in) | V_i, |
complex(dp), dimension(:,:,:), intent(out) | V_o, | ||
real(dp), dimension(:), intent(in) | r_i, | ||
real(dp), dimension(:), intent(in) | r_o, | ||
logical, intent(in) | extrap, | ||
integer, intent(out), optional | ivs_stat | ||
) |
Interpolation for a quantity V using splines.
Optionally, a variable 'ivs_stat' is returned that indicates whether the interpolation was a plain copy (1), a linear interpolation (2) or a spline interpolation (3).
Definition at line 801 of file sol_utilities.f90.
logical, public sol_utilities::debug_calc_xuq_arr = .false. |
plot debug information for calc_XUQ_arr
Definition at line 25 of file sol_utilities.f90.
logical, public sol_utilities::debug_interp_v_spline = .false. |
debug information for interp_v_spline
Definition at line 26 of file sol_utilities.f90.