PB3D
[2.45]
Ideal linear high-n MHD stability in 3-D
|
Numerical utilities related to PB3D operations. More...
Interfaces and Types | |
interface | conv_1d2nd |
Converts 1-D to n-D variables. More... | |
Functions/Subroutines | |
integer function, dimension(3), public | setup_par_id (grid, rich_lvl_max, rich_lvl_loc, tot_rich, par_lim, par_id_mem) |
Setup parallel id. More... | |
integer function, dimension(2), public | setup_rich_id (rich_lvl_max, tot_rich) |
Returns richardson id. More... | |
Numerical utilities related to PB3D operations.
integer function, dimension(3), public pb3d_utilities::setup_par_id | ( | type(grid_type), intent(in) | grid, |
integer, intent(in) | rich_lvl_max, | ||
integer, intent(in) | rich_lvl_loc, | ||
logical, intent(in), optional | tot_rich, | ||
integer, dimension(2), intent(in), optional | par_lim, | ||
integer, dimension(2), intent(inout), optional | par_id_mem | ||
) |
Setup parallel id.
The parallel id consists of:
par_id(1)
: start indexpar_id(2)
: end indexpar_id(3)
: strideThese are set up by considering the transformation between a point r \in (1\ldots n) in the local variable, with corresponding indices (a\ldots b) indicated by par_lim
in the HDF5 variable in memory, so that n = b-a+1:
p + k s = a + r - 1,
where k is an integer, s is the stride for Richardson level i (with max. level I):
\begin{aligned} s &= 2^{I-1} \ \text{for} \ i = 1 \\ &= 2^{I-i+1} \ \text{for} \ i > 1 \end{aligned}
and where p is the starting point in the HDF5 variables:
\begin{aligned} p &= 1 \ \text{for} \ i = 1 \\ &= 1 + s/2 \ \text{for} \ i > 1 \end{aligned}
Therefore, the minimum possible r lies in (0..s-1):
\mod(r-1,s) = r_\text{min}-1,
which leads to
r_\text{min} = 1 + \mod(p-a,s).
The maximum possible r, then, lies in (n-s+1..n), which leads to:
r_\text{max} = n - s + 1 + \mod(p-b-1,s).
These limits and the strides are saved in par_id
= \vec{r} = \begin{pmatrix}r_\text{min}\\ r_\text{max}\end{pmatrix}.
If the optional indices a and b are not given they are assumed to be 1 and n, with n = 1+ks, which simplifies the equations to:
\begin{aligned} r_\text{min} &= p \\ r_\text{max} &= n-p+1 . \end{aligned}
Optionally, the indices in the HDF5 arrays can also be returned in par_id_mem
= \begin{pmatrix}R_\text{min}\\ R_\text{max}\end{pmatrix}. These are equal to k-1, where k is the integer refered to above:
\vec{R} = 1 + \frac{a-1-p+\vec{r}}{s},
where the addition between a vector and a scalar should be seen as the element-wise operation.
[in] | rich_lvl_max | maximum Richardson level |
[in] | rich_lvl_loc | local Richardson level |
[in] | tot_rich | whether to combine with previous Richardson levels |
[in] | par_lim | limits (a..b) of variable requested |
[in,out] | par_id_mem | parallel id in memory |
Definition at line 181 of file PB3D_utilities.f90.
integer function, dimension(2), public pb3d_utilities::setup_rich_id | ( | integer, intent(in) | rich_lvl_max, |
logical, intent(in), optional | tot_rich | ||
) |
Returns richardson id.
rich_id
has the following structure:
rich_id(1)
: start Richardson levelrich_id(2)
: end Richardson level [in] | rich_lvl_max | maximum Richardson level |
[in] | tot_rich | whether to combine with previous Richardson levels |
Definition at line 234 of file PB3D_utilities.f90.