Extension: Geometries
#include "dg/geometries/geometries.h"
|
\( \nabla_\parallel f\) More...
Classes | |
struct | dg::geo::DS< ProductGeometry, IMatrix, Matrix, container > |
Class for the evaluation of parallel derivatives. More... | |
struct | dg::geo::Fieldaligned< ProductGeometry, IMatrix, container > |
Create and manage interpolation matrices from fieldline integration. More... | |
Typedefs | |
typedef ONE | dg::geo::FullLimiter |
Full Limiter means there is a limiter everywhere. More... | |
typedef ZERO | dg::geo::NoLimiter |
No Limiter. More... | |
Functions | |
template<class FieldAligned , class container > | |
void | dg::geo::assign_bc_along_field_2nd (const FieldAligned &fa, const container &fm, const container &f, const container &fp, container &fmg, container &fpg, dg::bc bound, std::array< double, 2 > boundary_value={0, 0}) |
Assign boundary conditions along magnetic field lines interpolating a 2nd order polynomial. More... | |
template<class FieldAligned , class container > | |
void | dg::geo::assign_bc_along_field_1st (const FieldAligned &fa, const container &fm, const container &fp, container &fmg, container &fpg, dg::bc bound, std::array< double, 2 > boundary_value={0, 0}) |
Assign boundary conditions along magnetic field lines interpolating a 1st order polynomial (a line) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::swap_bc_perp (const FieldAligned &fa, const container &fm, const container &fp, container &fmg, container &fpg) |
Swap the perp boundary condition. More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_forward (const FieldAligned &fa, double alpha, const container &f, const container &fp, double beta, container &g) |
forward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_forward2 (const FieldAligned &fa, double alpha, const container &f, const container &fp, const container &fpp, double beta, container &g) |
2nd order forward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_backward (const FieldAligned &fa, double alpha, const container &fm, const container &f, double beta, container &g) |
backward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_backward2 (const FieldAligned &fa, double alpha, const container &fmm, const container &fm, const container &f, double beta, container &g) |
2nd order backward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_centered (const FieldAligned &fa, double alpha, const container &fm, const container &fp, double beta, container &g) |
centered derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::dss_centered (const FieldAligned &fa, double alpha, const container &fm, const container &f, const container &fp, double beta, container &g) |
Centered derivative \( g = \alpha (\vec v\cdot \nabla)^2 f + \beta g \). More... | |
template<class FieldAligned , class container > | |
void | dg::geo::dssd_centered (const FieldAligned &fa, double alpha, const container &fm, const container &f, const container &fp, double beta, container &g) |
Centered derivative \( g = \alpha \nabla\cdot(\vec v \vec v\cdot \nabla) f + \beta g \). More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_divBackward (const FieldAligned &fa, double alpha, const container &fm, const container &f, double beta, container &g) |
backward derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_divForward (const FieldAligned &fa, double alpha, const container &f, const container &fp, double beta, container &g) |
forward derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_divCentered (const FieldAligned &fa, double alpha, const container &fm, const container &fp, double beta, container &g) |
centered derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\) More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_average (const FieldAligned &fa, double alpha, const container &fm, const container &fp, double beta, container &g) |
Compute average along a fieldline \( g = \alpha \frac{f_{k+1} + f_{k-1}}{2} + \beta g\). More... | |
template<class FieldAligned , class container > | |
void | dg::geo::ds_slope (const FieldAligned &fa, double alpha, const container &fm, const container &fp, double beta, container &g) |
Compute simple slope along a fieldline \( g = \alpha v^\varphi\frac{f_{k+1} - f_{k-1}}{2\Delta\varphi} + \beta g\). More... | |
template<class BinaryOp , class UnaryOp > | |
thrust::host_vector< double > | dg::geo::fieldaligned_evaluate (const aProductGeometry3d &grid, const CylindricalVectorLvl0 &vec, const BinaryOp &binary, const UnaryOp &unary, unsigned p0, unsigned rounds, double eps=1e-5) |
Evaluate a 2d functor and transform to all planes along the fieldlines More... | |
template<class BinaryOp , class UnaryOp > | |
MPI_Vector< thrust::host_vector< double > > | dg::geo::fieldaligned_evaluate (const aProductMPIGeometry3d &grid, const CylindricalVectorLvl0 &vec, const BinaryOp &binary, const UnaryOp &unary, unsigned p0, unsigned rounds, double eps=1e-5) |
Evaluate a 2d functor and transform to all planes along the fieldlines (MPI Version) More... | |
\( \nabla_\parallel f\)
typedef ONE dg::geo::FullLimiter |
Full Limiter means there is a limiter everywhere.
typedef ZERO dg::geo::NoLimiter |
No Limiter.
enum dg::geo::whichMatrix |
Enum for the use in Fieldaligned.
enum dg::geo::whichMatrix |
Enum for the use in Fieldaligned.
enum dg::geo::whichMatrix |
Enum for the use in Fieldaligned.
void dg::geo::assign_bc_along_field_1st | ( | const FieldAligned & | fa, |
const container & | fm, | ||
const container & | fp, | ||
container & | fmg, | ||
container & | fpg, | ||
dg::bc | bound, | ||
std::array< double, 2 > | boundary_value = {0,0} |
||
) |
Assign boundary conditions along magnetic field lines interpolating a 1st order polynomial (a line)
FieldAligned | |
container |
fa | |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
fmg | result (can alias fm) |
fpg | result (can alias fp) |
bound | either dg::NEU or dg::DIR (rest not implemented yet) |
boundary_value | first value is for incoming fieldlines, second one for outgoing |
dg::geo::Fieldaligned
object and the ones implicit in the einsPlus and einsMinus interpolations (Independently of whether the magnetic field was periodified or not). In this way it is possible to apply the same Fieldaligned object to several quantities with different boundary conditions and save quite a bit of memory consumption. void dg::geo::assign_bc_along_field_2nd | ( | const FieldAligned & | fa, |
const container & | fm, | ||
const container & | f, | ||
const container & | fp, | ||
container & | fmg, | ||
container & | fpg, | ||
dg::bc | bound, | ||
std::array< double, 2 > | boundary_value = {0,0} |
||
) |
Assign boundary conditions along magnetic field lines interpolating a 2nd order polynomial.
Call this function before one of the freestanding ds functions to replace the default boundary conditions with along field boundary conditions
This function replaces the values of the plus and minus fields that are outside of the domain with ghost values. These values are constructed by fitting a polynomial through the boundary point, the plus or minus point and the center point. For the exact resulting formula consult the parallel derivative writeup. This is achieved using masks that mark the points where fieldlines intersect the domain boundary and replace the interpolated boundary values.
FieldAligned | |
container |
fa | this object will be used to get grid distances and masking regions |
fm | fieldaligned(einsMinus, f, fm) |
f | |
fp | fieldaligned(einsPlus, f, fp) |
fmg | resulting eMinus field (can alias fm) |
fpg | resulting ePlus field (can alias fp) |
bound | either dg::NEU or dg::DIR (rest not implemented yet) |
boundary_value | first value is for incoming fieldlines, second one for outgoing |
dg::geo::Fieldaligned
object and the ones implicit in the einsPlus and einsMinus interpolations (Independently of whether the magnetic field was periodified or not). In this way it is possible to apply the same Fieldaligned object to several quantities with different boundary conditions and save quite a bit of memory consumption. void dg::geo::ds_average | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
Compute average along a fieldline \( g = \alpha \frac{f_{k+1} + f_{k-1}}{2} + \beta g\).
fa | this object will be used to get grid distances |
alpha | Scalar |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_backward | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | f, | ||
double | beta, | ||
container & | g | ||
) |
backward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\)
backward derivative \( g_i = \alpha \frac{v^\varphi}{\Delta\varphi}(f_{i} - f_{i-1}) + \beta g_i \)
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fm | fieldaligned(einsMinus, f, fm) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_backward2 | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fmm, | ||
const container & | fm, | ||
const container & | f, | ||
double | beta, | ||
container & | g | ||
) |
2nd order backward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\)
backward derivative \( g_i = \alpha \frac{v^\varphi}{2\Delta\varphi}(3f_{i} - 4f_{i-1} + f_{i-2}) + \beta g_i \)
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fm | fieldaligned(einsMinus, f, fm) |
fmm | twice apply fieldaligned(einsMinus, f, fm) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_centered | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
centered derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\)
The formula used is \( g_i = \alpha \frac{v^\varphi}{2\Delta\varphi}(f_{i+1} - f_{i-1}) + \beta g_i \)
fa | this object will be used to get grid distances |
alpha | Scalar |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_divBackward | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | f, | ||
double | beta, | ||
container & | g | ||
) |
backward derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\)
backward derivative \( g_i = \alpha \frac{1}{\Delta\varphi\sqrt{G_i}}(\sqrt{G_{i}}v^\varphi_{i}f_{i} - \sqrt{G_{i-1}}v^\varphi_{i-1}f_{i-1}) + \beta g_i\)
fa | this object will be used to get grid distances |
alpha | Scalar |
fm | fieldaligned(einsMinus, f, fm) |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_divCentered | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
centered derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\)
centered derivative \( g_i = \alpha \frac{1}{2\Delta\varphi\sqrt{G_i}}(\sqrt{G_{i+1}}v^\varphi_{i+1}f_{i+1} - \sqrt{G_{i-1}}v^\varphi_{i-1}f_{i-1}) + \beta g_i\)
fa | this object will be used to get grid distances |
alpha | Scalar |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_divForward | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | f, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
forward derivative \( g = \alpha \nabla \cdot \vec v f + \beta g\)
forward derivative \( g_i = \alpha \frac{1}{\Delta\varphi\sqrt{G_i}}(\sqrt{G_{i+1}}v^\varphi_{i+1}f_{i+1} - \sqrt{G_{i}}v^\varphi_{i}f_{i}) + \beta g_i\)
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
fp | fieldaligned(einsPlus, f, fp) |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_forward | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | f, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
forward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\)
forward derivative \( g_i = \alpha \frac{v^\varphi}{\Delta\varphi}(f_{i+1} - f_{i}) + \beta g_i\)
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fp | fieldaligned(einsPlus, f, fp) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_forward2 | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | f, | ||
const container & | fp, | ||
const container & | fpp, | ||
double | beta, | ||
container & | g | ||
) |
2nd order forward derivative \( g = \alpha \vec v \cdot \nabla f + \beta g\)
forward derivative \( g_i = \alpha \frac{v^\varphi}{2\Delta\varphi}(-f_{i+2} + 4f_{i+1} - 3f_{i}) + \beta g_i\)
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fp | fieldaligned(einsPlus, f, fp) |
fpp | twice apply fieldaligned(einsPlus, f, fp) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::ds_slope | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
Compute simple slope along a fieldline \( g = \alpha v^\varphi\frac{f_{k+1} - f_{k-1}}{2\Delta\varphi} + \beta g\).
fa | this object will be used to get grid distances |
alpha | Scalar |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
beta | Scalar |
g | contains result on output (may alias input vectors) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp ds_centered
void dg::geo::dss_centered | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | f, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
Centered derivative \( g = \alpha (\vec v\cdot \nabla)^2 f + \beta g \).
The formula used is
\[ \nabla_\parallel^2 f = 2\left(\frac{f^+}{h_z^+ h_z^0} - \frac{f^0}{h_z^- h_z^+} + \frac{f^-}{h_z^-h_z^0}\right) \]
which is the second derivative of a 2nd order polynomial fitted through the plus, minus and centre points the boundary conditions are implemented by mirroring points perpendicular to the boundary, which has some drawbacks as to the numerical stability and toroidal resolution.
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp void dg::geo::dssd_centered | ( | const FieldAligned & | fa, |
double | alpha, | ||
const container & | fm, | ||
const container & | f, | ||
const container & | fp, | ||
double | beta, | ||
container & | g | ||
) |
Centered derivative \( g = \alpha \nabla\cdot(\vec v \vec v\cdot \nabla) f + \beta g \).
The formula used is
\[ \Delta_\parallel f = \nabla\cdot (\vec v \nabla_\parallel f )\]
fa | this object will be used to get grid distances |
alpha | Scalar |
f | The vector to derive |
beta | Scalar |
g | contains result on output (may alias input vectors) |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
dg::geo::DS
but you have to compute the einsPlus and einsMinus interpolations from dg::geo::Fieldaligned
yourself. The reasoning for this function is that you can re-use the latter interpolations if you compute for example both first and second derivative of a function or even if you compute it for difference boundary conditions using dg::geo::assign_bc_along_field_2nd, dg::geo::assign_bc_along_field_1st or dg::geo::swap_bc_perp thrust::host_vector< double > dg::geo::fieldaligned_evaluate | ( | const aProductGeometry3d & | grid, |
const CylindricalVectorLvl0 & | vec, | ||
const BinaryOp & | binary, | ||
const UnaryOp & | unary, | ||
unsigned | p0, | ||
unsigned | rounds, | ||
double | eps = 1e-5 |
||
) |
Evaluate a 2d functor and transform to all planes along the fieldlines
The algorithm does the equivalent of the following:
BinaryOp
on a 2d planeu
is the given UnarayOp
, i
in [0..r] is the round index and j
in [0..Nz] is the plane index and \(\Delta\varphi\) is the angular distance given in the constructor (can be different from the actual grid distance hz!).j
, where the minus transformations get the inverted index \( N_z - j\).BinaryOp | Binary Functor |
UnaryOp | Unary Functor |
grid | The grid on which to integrate fieldlines. |
vec | The vector field to integrate. Note that you can control how the boundary conditions are represented by changing vec outside the grid domain using e.g. the periodify function. |
binary | Functor to evaluate in x-y |
unary | Functor to evaluate in z |
unary
is evaluated such that p0
corresponds to z=0, p0+1 corresponds to z=hz, p0-1 to z=-hz, ... p0 | The index of the plane to start |
rounds | The number of rounds r to follow a fieldline; can be zero, then the fieldlines are only followed within the current box ( no periodicity) |
eps | The accuracy of the fieldline integrator |
MPI_Vector< thrust::host_vector< double > > dg::geo::fieldaligned_evaluate | ( | const aProductMPIGeometry3d & | grid, |
const CylindricalVectorLvl0 & | vec, | ||
const BinaryOp & | binary, | ||
const UnaryOp & | unary, | ||
unsigned | p0, | ||
unsigned | rounds, | ||
double | eps = 1e-5 |
||
) |
Evaluate a 2d functor and transform to all planes along the fieldlines (MPI Version)
The algorithm does the equivalent of the following:
BinaryOp
on a 2d planeu
is the given UnarayOp
, i
in [0..r] is the round index and j
in [0..Nz] is the plane index and \(\Delta\varphi\) is the angular distance given in the constructor (can be different from the actual grid distance hz!).j
, where the minus transformations get the inverted index \( N_z - j\).BinaryOp | Binary Functor |
UnaryOp | Unary Functor |
grid | The grid on which to integrate fieldlines. |
vec | The vector field to integrate. Note that you can control how the boundary conditions are represented by changing vec outside the grid domain using e.g. the periodify function. |
binary | Functor to evaluate in x-y |
unary | Functor to evaluate in z |
unary
is evaluated such that p0
corresponds to z=0, p0+1 corresponds to z=hz, p0-1 to z=-hz, ... p0 | The index of the plane to start |
rounds | The number of rounds r to follow a fieldline; can be zero, then the fieldlines are only followed within the current box ( no periodicity) |
eps | The accuracy of the fieldline integrator |
void dg::geo::swap_bc_perp | ( | const FieldAligned & | fa, |
const container & | fm, | ||
const container & | fp, | ||
container & | fmg, | ||
container & | fpg | ||
) |
Swap the perp boundary condition.
This function multiplies (-1) to every value that lies outside the box. This effectively swaps the boundary conditions in the Fourier boundary mode, i.e. if NEU was used in Fieldaligned, then now they are DIR and vice versa.
FieldAligned | |
container |
fa | this object will be used to get masking regions |
fm | fieldaligned(einsMinus, f, fm) |
fp | fieldaligned(einsPlus, f, fp) |
fmg | resulting eMinus field (can alias fm) |
fpg | resulting ePlus field (can alias fp) |