PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Operations that have to do with the grids and different coordinate systems. More...
Functions/Subroutines | |
integer function, public | calc_norm_range (style, in_limits, eq_limits, X_limits, sol_limits, r_F_eq, r_F_X, r_F_sol, jq) |
Calculates normal range for the input grid, the equilibrium grid and/or the solution grid. More... | |
integer function, public | setup_grid_eq (grid_eq, eq_limits) |
Sets up the equilibrium grid. More... | |
integer function, public | setup_grid_eq_b (grid_eq, grid_eq_B, eq, only_half_grid) |
Sets up the field-aligned equilibrium grid. More... | |
integer function, public | setup_grid_x (grid_eq, grid_X, r_F_X, X_limits) |
Sets up the general perturbation grid, in which the perturbation variables are calculated. More... | |
integer function, public | setup_grid_sol (grid_eq, grid_X, grid_sol, r_F_sol, sol_limits) |
Sets up the general solution grid, in which the solution variables are calculated. More... | |
integer function, public | calc_ang_grid_eq_b (grid_eq, eq, only_half_grid) |
Calculate equilibrium grid that follows magnetic field lines. More... | |
integer function, public | redistribute_output_grid (grid, grid_out, no_outer_trim) |
Redistribute a grid to match the normal distribution of solution grid. More... | |
integer function, public | magn_grid_plot (grid) |
Plots the grid in real 3-D space. More... | |
integer function, public | print_output_grid (grid, grid_name, data_name, rich_lvl, par_div, remove_previous_arrs) |
Print grid variables to an output file. More... | |
Variables | |
logical, public | debug_calc_ang_grid_eq_b = .false. |
plot debug information for calc_ang_grid_eq_b() More... | |
Operations that have to do with the grids and different coordinate systems.
integer function, public grid_ops::calc_ang_grid_eq_b | ( | type(grid_type), intent(inout) | grid_eq, |
type(eq_1_type), intent(in), target | eq, | ||
logical, intent(in), optional | only_half_grid | ||
) |
Calculate equilibrium grid that follows magnetic field lines.
This grid is different from the equilibrium grid from setup_grid_eq for HELENA, as the latter is the output grid from HELENA, which is situated in a single poloidal cross-section, as opposed to a really field-aligne grid.
For VMEC, this is not used as the grid is field-aligned from the start.
only_half_grid
, only the even points of the parallel grid are calculated, which is useful for higher Richardson levels with VMEC so that only new angular points are calculated and the old ones reused.[in,out] | grid_eq | equilibrium grid of which to calculate angular part |
[in] | eq | flux equilibrium variables |
[in] | only_half_grid | calculate only half grid with even points |
Definition at line 798 of file grid_ops.f90.
integer function, public grid_ops::calc_norm_range | ( | character(len=*), intent(in) | style, |
integer, dimension(2), intent(inout), optional | in_limits, | ||
integer, dimension(2), intent(inout), optional | eq_limits, | ||
integer, dimension(2), intent(inout), optional | X_limits, | ||
integer, dimension(2), intent(inout), optional | sol_limits, | ||
real(dp), dimension(:), intent(inout), optional | r_F_eq, | ||
real(dp), dimension(:), intent(inout), optional, allocatable | r_F_X, | ||
real(dp), dimension(:), intent(inout), optional | r_F_sol, | ||
real(dp), dimension(:), intent(in), optional | jq | ||
) |
Calculates normal range for the input grid, the equilibrium grid and/or the solution grid.
General workings, depending on X_grid_style:
1 (equilibrium) | 2 (solution ) | 3 (enriched) |
---|---|---|
calc_eq() | calc_eq() | calc_eq() |
redistribute to sol | add points to optimize | |
interpolate to sol | interpolate to X | |
calc_x() | calc_x() | calc_x() |
redistribute to sol | copy to sol | redistribute to sol |
interpolate to sol | interpolate to sol |
[in] | style | style of calculation (PB3D: in, eq, X or sol; POST) |
[in,out] | in_limits | min. and max. index of in grid |
[in,out] | eq_limits | min. and max. index of eq grid for this process |
[in,out] | x_limits | min. and max. index of X grid for this process |
[in,out] | sol_limits | min. and max. index of sol grid for this process |
[in,out] | r_f_eq | equilibrium r_F |
[in,out] | r_f_x | perturbation r_F |
[in,out] | r_f_sol | solution r_F |
[in] | jq | q_saf (pol. flux) or rot_t (tor. flux) in total Flux variables |
Definition at line 46 of file grid_ops.f90.
integer function, public grid_ops::magn_grid_plot | ( | type(grid_type), intent(in) | grid | ) |
Plots the grid in real 3-D space.
This creates an animation that can be used by ParaView or VisIt.
The equilibrium grid should contain the fieldline-oriented angles with ang_1
the parallel angle and ang_2
the field line label.
ang_1
and ang_2
.n_theta_plot
and n_zeta_plot
from num_vars, but instead temporarily overwrites them with its own, since it is suposed to be 3-D also in the axisymmetric case.[in] | grid | fieldline-oriented equilibrium grid |
Definition at line 1115 of file grid_ops.f90.
integer function, public grid_ops::print_output_grid | ( | type(grid_type), intent(in) | grid, |
character(len=*), intent(in) | grid_name, | ||
character(len=*), intent(in) | data_name, | ||
integer, intent(in), optional | rich_lvl, | ||
logical, intent(in), optional | par_div, | ||
logical, intent(in), optional | remove_previous_arrs | ||
) |
Print grid variables to an output file.
If rich_lvl
is provided, _R_[rich_lvl]
is appended to the data name if it is > 0.
Optionally, it can be specified that this is a divided parallel grid, corresponding to the variable eq_jobs_lims
with index eq_job_nr
. In this case, the total grid size is adjusted to the one specified by eq_jobs_lims
and the grid is written as a subset.
grid_
is added in front the data_name.[in] | grid | grid variables |
[in] | grid_name | name to display |
[in] | data_name | name under which to store |
[in] | rich_lvl | Richardson level to reconstruct |
[in] | par_div | is a parallely divided grid |
[in] | remove_previous_arrs | remove previous variables if present |
Definition at line 1489 of file grid_ops.f90.
integer function, public grid_ops::redistribute_output_grid | ( | type(grid_type), intent(in) | grid, |
type(grid_type), intent(inout) | grid_out, | ||
logical, intent(in), optional | no_outer_trim | ||
) |
Redistribute a grid to match the normal distribution of solution grid.
The routine first calculates the smallest eq range that comprises the sol range. Then, it gets the lowest equilibrium limits able to setup an output grid that starts at index 1. After determining the output grid, it then sends the variables to their new processes using MPI.
no_outer_trim
.[in] | grid | equilibrium grid variables |
[in,out] | grid_out | redistributed equilibrium grid variables |
[in] | no_outer_trim | do not trim the outer limits |
Definition at line 995 of file grid_ops.f90.
integer function, public grid_ops::setup_grid_eq | ( | type(grid_type), intent(inout) | grid_eq, |
integer, dimension(2), intent(in) | eq_limits | ||
) |
Sets up the equilibrium grid.
The following variables are calculated:
For the implementation of the equilibrium grid the normal part of the grid is always given by the output of the equilibrium code, but the angular part depends on the style:
[in,out] | grid_eq | equilibrium grid |
[in] | eq_limits | min. and max. index of eq grid of this process |
Definition at line 529 of file grid_ops.f90.
integer function, public grid_ops::setup_grid_eq_b | ( | type(grid_type), intent(inout) | grid_eq, |
type(grid_type), intent(inout) | grid_eq_B, | ||
type(eq_1_type), intent(in) | eq, | ||
logical, intent(in), optional | only_half_grid | ||
) |
Sets up the field-aligned equilibrium grid.
This serves as a bridge to the solution grid, as it contains the same normal coordinate as the general grid, but the angular coordinates are defined by the solution grid.
Optionally, only half the grid can be calculated (i.e. only the even points), which is used for Richardson levels greater than 1.
setup_grid_eq
, the angular coordinates are also calculated here.[in,out] | grid_eq | general equilibrium grid |
[in,out] | grid_eq_b | field-aligned equilibrium grid |
[in] | eq | flux equilibrium variables |
[in] | only_half_grid | calculate only half grid with even points |
Definition at line 601 of file grid_ops.f90.
integer function, public grid_ops::setup_grid_sol | ( | type(grid_type), intent(in) | grid_eq, |
type(grid_type), intent(in) | grid_X, | ||
type(grid_type), intent(inout) | grid_sol, | ||
real(dp), dimension(:), intent(in) | r_F_sol, | ||
integer, dimension(2), intent(in) | sol_limits | ||
) |
Sets up the general solution grid, in which the solution variables are calculated.
For the solution grid, only one parallel point is used, but possibly multiple geodesic points, equal to the number of field lines, n_alpha
.
[in] | grid_eq | equilibrium grid |
[in] | grid_x | perturbation grid |
[in,out] | grid_sol | solution grid |
[in] | r_f_sol | points of solution grid |
[in] | sol_limits | min. and max. index of sol grid of this process |
Definition at line 728 of file grid_ops.f90.
integer function, public grid_ops::setup_grid_x | ( | type(grid_type), intent(in) | grid_eq, |
type(grid_type), intent(inout) | grid_X, | ||
real(dp), dimension(:), intent(in), allocatable | r_F_X, | ||
integer, dimension(2), intent(in) | X_limits | ||
) |
Sets up the general perturbation grid, in which the perturbation variables are calculated.
For X_grid_style
1, this grid is identical to the equilibrium grid, and for X_grid_style
2,it has the same angular extent but with different normal points, indicated by the variable r_F_X
.
[in] | grid_eq | equilibrium grid |
[in,out] | grid_x | perturbation grid |
[in] | r_f_x | points of perturbation grid |
[in] | x_limits | min. and max. index of perturbation grid of this process |
Definition at line 655 of file grid_ops.f90.
logical, public grid_ops::debug_calc_ang_grid_eq_b = .false. |
plot debug information for calc_ang_grid_eq_b()
Definition at line 25 of file grid_ops.f90.