WarpX
Loading...
Searching...
No Matches
CartesianNodalAlgorithm Struct Reference

#include <CartesianNodalAlgorithm.H>

Static Public Member Functions

static void InitializeStencilCoefficients (std::array< amrex::Real, 3 > &cell_size, amrex::Vector< amrex::Real > &stencil_coefs_x, amrex::Vector< amrex::Real > &stencil_coefs_y, amrex::Vector< amrex::Real > &stencil_coefs_z)
 
static amrex::Real ComputeMaxDt (amrex::Real const *const dx)
 
static amrex::IntVect GetMaxGuardCell ()
 Returns maximum number of guard cells required by the field-solve.
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDx (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDx (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_x, int const n_coefs_x, int const i, int const j, int const k, int const ncomp=0)
 
template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dxx (T_Field const &F, amrex::Real const *const coefs_x, int const, int const i, int const j, int const k, int const ncomp=0)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDy (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_y, int const, int const i, int const j, int const k, int const ncomp=0)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDy (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_y, int const n_coefs_y, int const i, int const j, int const k, int const ncomp=0)
 
template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dyy (T_Field const &F, amrex::Real const *const coefs_y, int const, int const i, int const j, int const k, int const ncomp=0)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real UpwardDz (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
 
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real DownwardDz (amrex::Array4< amrex::Real const > const &F, amrex::Real const *const coefs_z, int const n_coefs_z, int const i, int const j, int const k, int const ncomp=0)
 
template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real Dzz (T_Field const &F, amrex::Real const *const coefs_z, int const, int const i, int const j, int const k, int const ncomp=0)
 

Detailed Description

This struct contains only static functions to initialize the stencil coefficients and to compute finite-difference derivatives for the Cartesian nodal algorithm.

Member Function Documentation

◆ ComputeMaxDt()

static amrex::Real CartesianNodalAlgorithm::ComputeMaxDt ( amrex::Real const *const dx)
inlinestatic

Compute the maximum timestep, for which the scheme remains stable (Courant-Friedrichs-Levy limit)

◆ DownwardDx()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::DownwardDx ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_x,
int const n_coefs_x,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along x (For a solver on a staggered grid, UpwardDx and DownwardDx take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)

◆ DownwardDy()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::DownwardDy ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_y,
int const n_coefs_y,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along y (For a solver on a staggered grid, UpwardDy and DownwardDy take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)

◆ DownwardDz()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::DownwardDz ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_z,
int const n_coefs_z,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along z (For a solver on a staggered grid, UpwardDz and DownwardDz take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)

◆ Dxx()

template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::Dxx ( T_Field const & F,
amrex::Real const *const coefs_x,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform second derivative along x on any grid, from the field F

◆ Dyy()

template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::Dyy ( T_Field const & F,
amrex::Real const *const coefs_y,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform second derivative along y on any grid, from the field F

◆ Dzz()

template<typename T_Field>
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::Dzz ( T_Field const & F,
amrex::Real const *const coefs_z,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform second derivative along z on any grid, from the field F

◆ GetMaxGuardCell()

static amrex::IntVect CartesianNodalAlgorithm::GetMaxGuardCell ( )
inlinestatic

Returns maximum number of guard cells required by the field-solve.

◆ InitializeStencilCoefficients()

static void CartesianNodalAlgorithm::InitializeStencilCoefficients ( std::array< amrex::Real, 3 > & cell_size,
amrex::Vector< amrex::Real > & stencil_coefs_x,
amrex::Vector< amrex::Real > & stencil_coefs_y,
amrex::Vector< amrex::Real > & stencil_coefs_z )
inlinestatic

◆ UpwardDx()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::UpwardDx ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_x,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along x (For a solver on a staggered grid, UpwardDx and DownwardDx take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)

◆ UpwardDy()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::UpwardDy ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_y,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along y (For a solver on a staggered grid, UpwardDy and DownwardDy take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)

◆ UpwardDz()

AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE amrex::Real CartesianNodalAlgorithm::UpwardDz ( amrex::Array4< amrex::Real const > const & F,
amrex::Real const *const coefs_z,
int const ,
int const i,
int const j,
int const k,
int const ncomp = 0 )
inlinestatic

Perform derivative along z (For a solver on a staggered grid, UpwardDz and DownwardDz take into account the staggering; but for CartesianNodalAlgorithm, they are equivalent)


The documentation for this struct was generated from the following file: