File coordinates.hxx#

class Coordinates#
#include <coordinates.hxx>

Represents a coordinate system, and associated operators

This is a container for a collection of metric tensor components

Public Types

using FieldMetric = Field3D#

Public Functions

Coordinates(Mesh *mesh, Options *options = nullptr)#

Standard constructor from input.

Coordinates(Mesh *mesh, Options *options, const CELL_LOC loc, const Coordinates *coords_in, bool force_interpolate_from_centre = false)#

Constructor interpolating from another Coordinates object By default attempts to read staggered Coordinates from grid data source, interpolating from CELL_CENTRE if not present. Set force_interpolate_from_centre argument to true to always interpolate (useful if CELL_CENTRE Coordinates have been changed, so reading from file would not be correct).

Coordinates(Mesh *mesh, FieldMetric dx, FieldMetric dy, FieldMetric dz, FieldMetric J, FieldMetric Bxy, FieldMetric g11, FieldMetric g22, FieldMetric g33, FieldMetric g12, FieldMetric g13, FieldMetric g23, FieldMetric g_11, FieldMetric g_22, FieldMetric g_33, FieldMetric g_12, FieldMetric g_13, FieldMetric g_23, FieldMetric ShiftTorsion, FieldMetric IntShiftTorsion)#

A constructor useful for testing purposes. To use it, inherit from Coordinates. If calculate_geometry is true (default), calculate the non-uniform variables, Christoffel symbols

Coordinates &operator=(Coordinates&&) = default#
~Coordinates() = default#
void outputVars(Options &output_options)#

Add variables to output_options, for post-processing.

const Field2D &zlength() const#

Length of the Z domain. Used for FFTs.

int geometry(bool recalculate_staggered = true, bool force_interpolate_from_centre = false)#

Calculate differential geometry quantities from the metric tensor.

int calcCovariant(const std::string &region = "RGN_ALL")#

Invert contravatiant metric to get covariant components.

int calcContravariant(const std::string &region = "RGN_ALL")#

Invert covariant metric to get contravariant components.

int jacobian()#

Calculate J and Bxy.

inline void setParallelTransform(std::unique_ptr<ParallelTransform> pt)#

Set the parallel (y) transform for this mesh. Mostly useful for tests.

inline ParallelTransform &getParallelTransform()#

Return the parallel transform.

FieldMetric DDX(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY")#
FieldMetric DDY(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY") const#
FieldMetric DDZ(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY")#
Field3D DDX(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY")#
Field3D DDY(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY") const#
Field3D DDZ(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT", const std::string &region = "RGN_NOBNDRY")#
FieldMetric Grad_par(const Field2D &var, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#

Gradient along magnetic field b.Grad(f)

Field3D Grad_par(const Field3D &var, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#
FieldMetric Vpar_Grad_par(const Field2D &v, const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#

Advection along magnetic field V*b.Grad(f)

Field3D Vpar_Grad_par(const Field3D &v, const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#
FieldMetric Div_par(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#

Divergence along magnetic field Div(b*f) = B.Grad(f/B)

Field3D Div_par(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#
FieldMetric Grad2_par2(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#
Field3D Grad2_par2(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &method = "DEFAULT")#
FieldMetric Delp2(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, bool useFFT = true)#
Field3D Delp2(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, bool useFFT = true)#
FieldPerp Delp2(const FieldPerp &f, CELL_LOC outloc = CELL_DEFAULT, bool useFFT = true)#
FieldMetric Laplace_par(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT)#
Field3D Laplace_par(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT)#
FieldMetric Laplace(const Field2D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &dfdy_boundary_conditions = "free_o3", const std::string &dfdy_dy_region = "")#
Field3D Laplace(const Field3D &f, CELL_LOC outloc = CELL_DEFAULT, const std::string &dfdy_boundary_conditions = "free_o3", const std::string &dfdy_dy_region = "")#
Field2D Laplace_perpXY(const Field2D &A, const Field2D &f)#

Public Members

FieldMetric dx#
FieldMetric dy#
FieldMetric dz#

Mesh spacing in x, y and z.

bool non_uniform#

True if corrections for non-uniform mesh spacing should be included in operators.

FieldMetric d1_dx#

2nd-order correction for non-uniform meshes d/di(1/dx), d/di(1/dy) and d/di(1/dz)

FieldMetric d1_dy#
FieldMetric d1_dz#
FieldMetric J#

Coordinate system Jacobian, so volume of cell is J*dx*dy*dz.

FieldMetric Bxy#

Magnitude of B = nabla z times nabla x.

FieldMetric g11#

Contravariant metric tensor (g^{ij})

FieldMetric g22#
FieldMetric g33#
FieldMetric g12#
FieldMetric g13#
FieldMetric g23#
FieldMetric g_11#

Covariant metric tensor.

FieldMetric g_22#
FieldMetric g_33#
FieldMetric g_12#
FieldMetric g_13#
FieldMetric g_23#
FieldMetric G1_11#

Christoffel symbol of the second kind (connection coefficients)

FieldMetric G1_22#
FieldMetric G1_33#
FieldMetric G1_12#
FieldMetric G1_13#
FieldMetric G1_23#
FieldMetric G2_11#
FieldMetric G2_22#
FieldMetric G2_33#
FieldMetric G2_12#
FieldMetric G2_13#
FieldMetric G2_23#
FieldMetric G3_11#
FieldMetric G3_22#
FieldMetric G3_33#
FieldMetric G3_12#
FieldMetric G3_13#
FieldMetric G3_23#
FieldMetric G1#
FieldMetric G2#
FieldMetric G3#
FieldMetric ShiftTorsion#

d pitch angle / dx. Needed for vector differentials (Curl)

FieldMetric IntShiftTorsion#

Integrated shear (I in BOUT notation)

Private Functions

void setParallelTransform(Options *options)#

Set the parallel (y) transform from the options file. Used in the constructor to create the transform object.

const FieldMetric &invSg() const#
const FieldMetric &Grad2_par2_DDY_invSg(CELL_LOC outloc, const std::string &method) const#
void checkCovariant()#
void checkContravariant()#

Private Members

int nz#
Mesh *localmesh#
CELL_LOC location#
std::unique_ptr<ParallelTransform> transform = {nullptr}#

Handles calculation of yup and ydown.

mutable std::unique_ptr<Field2D> zlength_cache = {nullptr}#

Cache variable for zlength. Invalidated when Coordinates::geometry is called

mutable std::map<std::string, std::unique_ptr<FieldMetric>> Grad2_par2_DDY_invSgCache#

Cache variable for Grad2_par2.

mutable std::unique_ptr<FieldMetric> invSgCache = {nullptr}#