PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Operations and variables pertaining to the vacuum response. More...
Functions/Subroutines | |
integer function, public | store_vac (grid, eq_1, eq_2, vac) |
Stores calculation of the vacuum response by storing the necessary variables. More... | |
integer function | calc_gh (vac, n_r_in, lims_r_in, x_vec_in, G_in, H_in) |
Calculate the matrices G and H . More... | |
integer function | calc_gh_1 (vac, n_r_in, lims_r_in, x_vec_in, G_in, H_in, ext_in) |
Calculates matrices G and H in 3-D configuration. More... | |
integer function | calc_gh_2 (vac, n_r_in, lims_r_in, x_vec_in, G_in, H_in, ext_in) |
Calculates matrices G and H in axisymmetric configurations. More... | |
integer function, public | calc_vac_res (mds, vac) |
Calculates the vacuum response. More... | |
integer function, public | vac_pot_plot (grid_sol, vac, sol, X_id) |
Calculate vacuum potential at some positions that are not resonant. More... | |
integer function | solve_phi_bem (vac, R, Phi, n_RPhi, n_loc_RPhi, lims_c_RPhi, desc_RPhi) |
Calculates potential Phi on boundary of BEM in terms of some right-hand side . More... | |
integer function, public | print_output_vac (vac, data_name, rich_lvl) |
Print vacuum quantities of a certain order to an output file. More... | |
Variables | |
logical, public | debug_calc_gh = .false. |
plot debug information for calc_GH() More... | |
logical, public | debug_calc_vac_res = .false. |
plot debug information for calc_vac_res() More... | |
logical, public | debug_vac_pot_plot = .false. |
plot debug information for vac_pot_plot() More... | |
Operations and variables pertaining to the vacuum response.
The vacuum response is calculated using the Boundary Element Method. There are two possibilities, indicated by the variable style
in vac_vars
.
The final matrix equation is solved through Strumpack [10].
integer function vac_ops::calc_gh | ( | type(vac_type), intent(inout), target | vac, |
integer, intent(in), optional | n_r_in, | ||
integer, dimension(:,:), intent(in), optional, target | lims_r_in, | ||
real(dp), dimension(:,:), intent(in), optional, target | x_vec_in, | ||
real(dp), dimension(:,:), intent(in), optional, target | G_in, | ||
real(dp), dimension(:,:), intent(in), optional, target | H_in | ||
) |
Calculate the matrices G
and H
.
This is done by iterating over the subintervals, calculating the contributions to the potential values on the left and right boundary of the interval, and adding these to the appropriate values in G and H.
As each process in the blacs grid only has access to its own values, an extra interval is is calculated to the left and right of the internal process grid edges.
Optionally, this procedure can be used to calculate the coefficients G and H with respect to other points, outside of the boundary. This is useful when the potential is to be calculated at external points.
In this case, however, no (near-) singular points are calculated, and if they appear any way, perhaps by accident, they will not be accurate.
\Note Multiple field lines are stored sequentially, which implies that the integral between the last point on a field line and the first point on the next field lines needs to be left out.
[in,out] | vac | vacuum variables |
[in] | n_r_in | total number of rows of external poins of influence |
[in] | lims_r_in | row limits of external points of influence |
[in] | x_vec_in | external position of influence |
[in] | g_in | external G |
[in] | h_in | external H |
Definition at line 501 of file vac_ops.f90.
integer function vac_ops::calc_gh_1 | ( | type(vac_type), intent(inout) | vac, |
integer, intent(in) | n_r_in, | ||
integer, dimension(:,:), intent(in) | lims_r_in, | ||
real(dp), dimension(:,:), intent(in) | x_vec_in, | ||
real(dp), dimension(:,:), intent(in), pointer | G_in, | ||
real(dp), dimension(:,:), intent(in), pointer | H_in, | ||
logical, intent(in) | ext_in | ||
) |
Calculates matrices G
and H
in 3-D configuration.
[in,out] | vac | vacuum variables |
[in] | n_r_in | total number of rows of external poins of influence |
[in] | lims_r_in | row limits of external points of influence |
[in] | x_vec_in | position of influence |
[in] | g_in | G at position of influence |
[in] | h_in | H at position of influence |
[in] | ext_in | position of influence is external |
Definition at line 779 of file vac_ops.f90.
integer function vac_ops::calc_gh_2 | ( | type(vac_type), intent(inout) | vac, |
integer, intent(in) | n_r_in, | ||
integer, dimension(:,:), intent(in) | lims_r_in, | ||
real(dp), dimension(:,:), intent(in) | x_vec_in, | ||
real(dp), dimension(:,:), intent(in), pointer | G_in, | ||
real(dp), dimension(:,:), intent(in), pointer | H_in, | ||
logical, intent(in) | ext_in | ||
) |
Calculates matrices G
and H
in axisymmetric configurations.
It makes use of toroidal functions.
[in,out] | vac | vacuum variables |
[in] | n_r_in | total number of rows of external poins of influence |
[in] | lims_r_in | row limits of external points of influence |
[in] | x_vec_in | position of influence |
[in] | g_in | G at position of influence |
[in] | h_in | H at position of influence |
[in] | ext_in | position of influence is external |
Definition at line 983 of file vac_ops.f90.
integer function, public vac_ops::calc_vac_res | ( | type(modes_type), intent(in) | mds, |
type(vac_type), intent(inout) | vac | ||
) |
Calculates the vacuum response.
First G
and H
are completed if this is not the first Richardson level.
Then the vacuum contribution is calculated using the two-fold strategy of first solving \(\overline{\text{H}}\overline{\Phi} = \overline{\text{G}}\overline{\text{E}} \overline{\text{P}}\), followed by left-multiplication of \(\overline{\Phi}\) by \(\overline{\text{P}}\overline{\text{E}}^\dagger\).
The non-square matrix \(\overline{\text{E}}\) contains the exponents for different mode numbers and different poloidal grid points, whereas \(\overline{\text{P}}\) are diagonal matrices that contain the factors \((nq-m)\).
In practice, the complex matrix \(\overline{\text{E}}\) is split in the two components of a real matrix twice the width.
For vacuum style 1, the poloidal grid points correspond to the paralell grid points and have to be provied by a grid
variable.
If jump_to_sol
is used for the current Richardson level, the vacuum quantities are not calculated, but just restored.
Currently, this procedure only works for vacuum style 2 (axisymmetric).
[in] | mds | general modes variables |
[in,out] | vac | vacuum variables |
Definition at line 1417 of file vac_ops.f90.
integer function, public vac_ops::print_output_vac | ( | type(vac_type), intent(in) | vac, |
character(len=*), intent(in) | data_name, | ||
integer, intent(in), optional | rich_lvl | ||
) |
Print vacuum quantities of a certain order to an output file.
misc_vac
norm
x_vec
vac_res
If rich_lvl
is provided, "_R_[rich_lvl]"
is appended to the data name if it is >0
.
H
and G
matrices because they are large and foreseen to be saved directly as Hierarchical matrices in the future. They therefore need to be regenerated if Richardson restart is done, or when after a jump to solution, another Richardson level is attempted. In calc_vac, it is checked when copying the vacuum variables from previous to current Richardson level, whether the G
and H
matrices are allocated and if not, they are calculated.[in] | vac | vacuum variables |
[in] | data_name | name under which to store |
[in] | rich_lvl | Richardson level to print |
Definition at line 2347 of file vac_ops.f90.
integer function vac_ops::solve_phi_bem | ( | type(vac_type), intent(in), target | vac, |
real(dp), dimension(:,:), intent(in) | R, | ||
real(dp), dimension(:,:), intent(inout), target | Phi, | ||
integer, dimension(2), intent(in) | n_RPhi, | ||
integer, dimension(2), intent(in) | n_loc_RPhi, | ||
integer, dimension(:,:), intent(in) | lims_c_RPhi, | ||
integer, dimension(blacsctxtsize), intent(in), target | desc_RPhi | ||
) |
Calculates potential Phi
on boundary of BEM in terms of some right-hand side .
This is done by solving with STRUMPack the system of equations \(\overline{\text{H}}\overline{\Phi} = \overline{\text{G}}\overline{\text{R}}\), where \(\overline{\text{R}}\) is the right-hand side, for example equal to \(\overline{\text{E}}\overline{\text{P}}\) to solve in function of the individual Fourier modes.
[in] | vac | vacuum variables |
[in,out] | phi | \(\overline{\Phi}\) |
[in] | r | \(\overline{\text{R}}\) |
[in] | n_rphi | n for \(\overline{\text{R}}\) and \(\overline{\Phi}\) |
[in] | n_loc_rphi | local n for \(\overline{\text{R}}\) and \(\overline{\Phi}\) |
[in] | lims_c_rphi | local limits for row and columns of \(\overline{\text{R}}\) and \(\overline{\Phi}\) |
[in] | desc_rphi | descriptor for \(\overline{\text{R}}\) and \(\overline{\Phi}\) |
Definition at line 2163 of file vac_ops.f90.
integer function, public vac_ops::store_vac | ( | type(grid_type), intent(in) | grid, |
type(eq_1_type), intent(in) | eq_1, | ||
type(eq_2_type), intent(in) | eq_2, | ||
type(vac_type), intent(inout) | vac | ||
) |
Stores calculation of the vacuum response by storing the necessary variables.
This is done by the last process, which is the one that contains the variables at the plasma edge, and then this is broadcasted to the other processes.
The workings of this routine are summarized in the following diagram for VMEC:
Before storing the new variables, the procedure checks the following things:
If there is a jump straight to the solution, however, the procedure returns after reconstructing the current level's variables.
If the previous G
and H
variables are empty, they are regenerated later, in calc_gh(). This indicates that a reconstruction happened, either because a Richardson restart was performed, or because a new Richardson level was started after a previous level in which it was jumped straight to the solution.
For HELENA, there is only 1 equilibrium job at the first Richardson level, and the vacuum has to be calculated only once then. If there is a jump to solution or if there is Richardson restart, it only needs to be reconstructed.
[in] | grid | equilibrium grid |
[in] | eq_1 | flux equilibrium variables |
[in] | eq_2 | metric equilibrium variables |
[in,out] | vac | vacuum variables |
Definition at line 113 of file vac_ops.f90.
integer function, public vac_ops::vac_pot_plot | ( | type(grid_type), intent(in) | grid_sol, |
type(vac_type), intent(inout) | vac, | ||
type(sol_type), intent(in) | sol, | ||
integer, intent(in) | X_id | ||
) |
Calculate vacuum potential at some positions that are not resonant.
This is done by using the relation between H and G on the plasma boundary, \(\overline{\text{H}}_\text{s} \overline{\Phi}_\text{s} = \overline{\text{G}}_\text{s} \overline{\text{D}\phi}_\text{s}\), where \(\overline{\text{D}\phi}_\text{s}\) is the normal derivative of the potential on the plasma edge. The matrices \(\overline{\text{H}}\) and \(\overline{\text{G}}\) have to be calculated in advance.
This relation is then used at different positions to calculate the potential at a different position through \(4 \pi \phi = \left[ - \overline{\text{H}} \overline{\text{I}} \overline{\text{H}}_\text{s}^{-1} \overline{\text{G}}_\text{s} + \overline{\text{G}} \overline{\text{I}} \right] \overline{\text{D}\phi}_\text{s}\), where \(\overline{\text{I}}\) is an integration rule.
This system of equations contains one row per point at which the potential is to be calculated.
Currently, this routine creates a 3-D regular grid that is equidistant in cylindrical coordinates, with n_vac_plot
values in R and Z, and n_zeta_plot
values in the cylindrical angle. The limits in these variables, furthermore, are given by min_Rvac_plot
, max_Rvac_plot
, min_Zvac_plot
, max_Zvac_plot
, min_zeta_plot
and max_zeta_plot
.
[in] | grid_sol | solution grid |
[in,out] | vac | vacuum variables |
[in] | sol | solution variables |
[in] | x_id | nr. of Eigenvalue |
Definition at line 1857 of file vac_ops.f90.
logical, public vac_ops::debug_calc_gh = .false. |
plot debug information for calc_GH()
Definition at line 46 of file vac_ops.f90.
logical, public vac_ops::debug_calc_vac_res = .false. |
plot debug information for calc_vac_res()
Definition at line 47 of file vac_ops.f90.
logical, public vac_ops::debug_vac_pot_plot = .false. |
plot debug information for vac_pot_plot()
Definition at line 48 of file vac_ops.f90.