PB3D  [2.45]
Ideal linear high-n MHD stability in 3-D
Functions/Subroutines | Variables
vac_ops Module Reference

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...
 

Detailed Description

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.

  1. 3-D field aligned:
    • Makes use of the collocation method for the grid points along the magnetic field lines.
    • The 3-D integrals are therefore integrated using Weyl's theorem [6], p. 5.
  2. axisymmetric:
    • Makes use of an analytical expression for the toroidal integration of the Green's functions, as shown in [3].
    • The integration in the poloidal angle is done using the collocation method.

The final matrix equation is solved through Strumpack [10].

See also
See [17].
Todo:
The vacuum part of PB3D is still under construction and not usable yet.

Function/Subroutine Documentation

◆ calc_gh()

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.

Todo:
For 3-D vacuums, step sizes are constant. Subsequent Richardson levels should therefore make use of the fact that the contribution to the points inherited from the previous levels can be just scaled by 1/2 and do not need to be recalculated. Currently, the copy is done correctly in interlaced_vac_copy(), but they are afterwards overwritten.
Returns
ierr
Parameters
[in,out]vacvacuum variables
[in]n_r_intotal number of rows of external poins of influence
[in]lims_r_inrow limits of external points of influence
[in]x_vec_inexternal position of influence
[in]g_inexternal G
[in]h_inexternal H

Definition at line 501 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ calc_gh_1()

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.

See also
calc_gh.
Parameters
[in,out]vacvacuum variables
[in]n_r_intotal number of rows of external poins of influence
[in]lims_r_inrow limits of external points of influence
[in]x_vec_inposition of influence
[in]g_inG at position of influence
[in]h_inH at position of influence
[in]ext_inposition of influence is external

Definition at line 779 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ calc_gh_2()

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.

See also
calc_gh.
Parameters
[in,out]vacvacuum variables
[in]n_r_intotal number of rows of external poins of influence
[in]lims_r_inrow limits of external points of influence
[in]x_vec_inposition of influence
[in]g_inG at position of influence
[in]h_inH at position of influence
[in]ext_inposition of influence is external

Definition at line 983 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ calc_vac_res()

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).

Returns
ierr
Parameters
[in]mdsgeneral modes variables
[in,out]vacvacuum variables

Definition at line 1417 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ print_output_vac()

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.

  • vac:
    • 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.

Note
  1. This routine should be called by a single process, in contrast to the other output routines such as eq_ops.print_output_eq(), print_output_sol(), ...
  2. This process should be the last one, as it will set the boundary contribution.
  3. This routine does not save the 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.
Returns
ierr
Parameters
[in]vacvacuum variables
[in]data_namename under which to store
[in]rich_lvlRichardson level to print

Definition at line 2347 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve_phi_bem()

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.

Parameters
[in]vacvacuum variables
[in,out]phi\(\overline{\Phi}\)
[in]r\(\overline{\text{R}}\)
[in]n_rphin for \(\overline{\text{R}}\) and \(\overline{\Phi}\)
[in]n_loc_rphilocal n for \(\overline{\text{R}}\) and \(\overline{\Phi}\)
[in]lims_c_rphilocal limits for row and columns of \(\overline{\text{R}}\) and \(\overline{\Phi}\)
[in]desc_rphidescriptor for \(\overline{\text{R}}\) and \(\overline{\Phi}\)

Definition at line 2163 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ store_vac()

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:

  • For the first equilibrium job, this routine copies the results from the previous Richardson level into the appropriate subranges of the new vacuum variables if no restart of Richardson levels was done.
  • If a restart was done and the level is greater than one, there is a reconstruction of the previous level's variables. But this will happen in calc_vac_res(), not in this procedure.

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.

Returns
ierr
Parameters
[in]gridequilibrium grid
[in]eq_1flux equilibrium variables
[in]eq_2metric equilibrium variables
[in,out]vacvacuum variables

Definition at line 113 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ vac_pot_plot()

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.

Parameters
[in]grid_solsolution grid
[in,out]vacvacuum variables
[in]solsolution variables
[in]x_idnr. of Eigenvalue

Definition at line 1857 of file vac_ops.f90.

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Variable Documentation

◆ debug_calc_gh

logical, public vac_ops::debug_calc_gh = .false.

plot debug information for calc_GH()

Note
Debug version only

Definition at line 46 of file vac_ops.f90.

◆ debug_calc_vac_res

logical, public vac_ops::debug_calc_vac_res = .false.

plot debug information for calc_vac_res()

Note
Debug version only

Definition at line 47 of file vac_ops.f90.

◆ debug_vac_pot_plot

logical, public vac_ops::debug_vac_pot_plot = .false.

plot debug information for vac_pot_plot()

Note
Debug version only

Definition at line 48 of file vac_ops.f90.