ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TerrainMetrics.H File Reference
#include <AMReX.H>
#include <AMReX_Geometry.H>
#include <AMReX_MultiFab.H>
#include <ERF_IndexDefines.H>
Include dependency graph for ERF_TerrainMetrics.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void init_default_zphys (int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
 
void init_zlevels (amrex::Vector< amrex::Vector< amrex::Real >> &zlevels_stag, amrex::Vector< amrex::Vector< amrex::Real >> &stretched_dz_h, amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real >> &stretched_dz_d, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::IntVect > const &ref_ratio, const amrex::Real grid_stretching_ratio, const amrex::Real zsurf, const amrex::Real dz0)
 
void make_terrain_fitted_coords (int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::Vector< amrex::Real > const &z_levels_h, amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
 
void init_which_terrain_grid (int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::Vector< amrex::Real > const &z_levels_h)
 
AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real Compute_h_zeta_AtCellCenter (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtCellCenter (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtCellCenter (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtIface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtJface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtKface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtKface (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterK (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterK (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterK (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterJ (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterJ (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterJ (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterI (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterI (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterI (const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtCellCenter (const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace (const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Zrel_AtCellCenter (const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real OmegaFromW (int &i, int &j, int &k, amrex::Real w, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega (int &i, int &j, int &k, amrex::Real omega, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_scalar_flux (const int &i, const int &j, const int &klo, const amrex::Real &flux, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &zphys_arr, const amrex::Array4< amrex::Real > &phi1_arr, const amrex::Array4< amrex::Real > &phi2_arr, const amrex::Array4< amrex::Real > &phi3_arr)
 
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_stress_tensor (const int &i, const int &j, const int &klo, const amrex::Real &flux, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &zphys_arr, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &w_arr, const amrex::Array4< amrex::Real > &tau11_arr, const amrex::Array4< amrex::Real > &tau22_arr, const amrex::Array4< amrex::Real > &tau33_arr, const amrex::Array4< amrex::Real > &tau12_arr, const amrex::Array4< amrex::Real > &tau21_arr, const amrex::Array4< amrex::Real > &tau13_arr, const amrex::Array4< amrex::Real > &tau31_arr, const amrex::Array4< amrex::Real > &tau23_arr, const amrex::Array4< amrex::Real > &tau32_arr)
 

Function Documentation

◆ Compute_h_eta_AtCellCenter()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtCellCenter ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
80 {
81  amrex::Real dyInv = cellSizeInv[1];
82  amrex::Real met_h_eta = 0.25 * dyInv *
83  ( z_nd(i,j+1,k) + z_nd(i+1,j+1,k) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1)
84  -z_nd(i,j ,k) - z_nd(i+1,j ,k) - z_nd(i,j ,k+1) - z_nd(i+1,j ,k+1) );
85  return met_h_eta;
86 }
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by ComputeStrain_T(), ComputeStressConsVisc_T(), ComputeStressVarVisc_T(), ERFPhysBCFunct_cons::impose_vertical_cons_bcs(), rotate_scalar_flux(), MOSTAverage::set_norm_indices_T(), and MOSTAverage::set_norm_positions_T().

Here is the caller graph for this function:

◆ Compute_h_eta_AtEdgeCenterI()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterI ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
341 {
342  amrex::Real dyInv = cellSizeInv[1];
343  amrex::Real met_h_eta = 0.25 * dyInv *
344  ( z_nd(i+1,j+1,k) + z_nd(i,j+1,k)
345  -z_nd(i+1,j-1,k) - z_nd(i,j-1,k) );
346  return met_h_eta;
347 }

Referenced by ComputeStrain_T(), ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_eta_AtEdgeCenterJ()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterJ ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
299 {
300  amrex::Real dyInv = cellSizeInv[1];
301  amrex::Real met_h_eta = dyInv * ( z_nd(i,j+1,k) - z_nd(i,j,k) );
302  return met_h_eta;
303 }

Referenced by ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_eta_AtEdgeCenterK()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterK ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
253 {
254  amrex::Real dyInv = cellSizeInv[1];
255  amrex::Real met_h_eta = 0.25 * dyInv *
256  ( z_nd(i,j+1,k) + z_nd(i,j+1,k+1)
257  -z_nd(i,j-1,k) - z_nd(i,j-1,k+1) );
258  return met_h_eta;
259 }

Referenced by ComputeStrain_T().

Here is the caller graph for this function:

◆ Compute_h_eta_AtIface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtIface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
126 {
127  amrex::Real met_h_eta = 0.5 * cellSizeInv[1] * ( (z_nd(i,j+1,k ) - z_nd(i,j,k ))
128  + (z_nd(i,j+1,k+1) - z_nd(i,j,k+1)) );
129  return met_h_eta;
130 }

Referenced by ERFPhysBCFunct_u::impose_vertical_xvel_bcs(), and rotate_stress_tensor().

Here is the caller graph for this function:

◆ Compute_h_eta_AtJface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
165 {
166  amrex::Real met_h_eta = 0.125 * cellSizeInv[1] *
167  ( (z_nd(i ,j+1,k ) - z_nd(i ,j-1,k )) + (z_nd(i ,j+1,k+1) - z_nd(i ,j-1,k+1)) +
168  (z_nd(i+1,j+1,k ) - z_nd(i+1,j-1,k )) + (z_nd(i+1,j+1,k+1) - z_nd(i+1,j-1,k+1)) );
169  return met_h_eta;
170 }

Referenced by DiffusionSrcForState_T(), erf_fast_rhs_MT(), erf_fast_rhs_T(), ERFPhysBCFunct_v::impose_vertical_yvel_bcs(), and MOSTAverage::set_rotated_fields().

Here is the caller graph for this function:

◆ Compute_h_eta_AtKface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtKface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
206 {
207  amrex::Real met_h_eta = 0.5 * cellSizeInv[1] * ( (z_nd(i ,j+1,k) - z_nd(i ,j,k))
208  + (z_nd(i+1,j+1,k) - z_nd(i+1,j,k)) );
209  return met_h_eta;
210 }

Referenced by DiffusionSrcForState_T().

Here is the caller graph for this function:

◆ Compute_h_xi_AtCellCenter()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtCellCenter ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
65 {
66  amrex::Real dxInv = cellSizeInv[0];
67  amrex::Real met_h_xi = 0.25 * dxInv *
68  ( z_nd(i+1,j,k) + z_nd(i+1,j+1,k) + z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k+1)
69  -z_nd(i ,j,k) - z_nd(i ,j+1,k) - z_nd(i ,j,k+1) - z_nd(i ,j+1,k+1) );
70  return met_h_xi;
71 }

Referenced by ComputeStrain_T(), ComputeStressConsVisc_T(), ComputeStressVarVisc_T(), ERFPhysBCFunct_cons::impose_vertical_cons_bcs(), rotate_scalar_flux(), MOSTAverage::set_norm_indices_T(), and MOSTAverage::set_norm_positions_T().

Here is the caller graph for this function:

◆ Compute_h_xi_AtEdgeCenterI()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterI ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
328 {
329  amrex::Real dxInv = cellSizeInv[0];
330  amrex::Real met_h_xi = dxInv * ( z_nd(i+1,j,k) - z_nd(i,j,k) );
331  return met_h_xi;
332 }

Referenced by ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_xi_AtEdgeCenterJ()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterJ ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
284 {
285  amrex::Real dxInv = cellSizeInv[0];
286  amrex::Real met_h_xi = 0.25 * dxInv *
287  ( z_nd(i+1,j+1,k) + z_nd(i+1,j ,k)
288  -z_nd(i-1,j+1,k) - z_nd(i-1,j ,k) );
289  return met_h_xi;
290 }

Referenced by ComputeStrain_T(), ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_xi_AtEdgeCenterK()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterK ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
238 {
239  amrex::Real dxInv = cellSizeInv[0];
240  amrex::Real met_h_xi = 0.25 * dxInv *
241  ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1)
242  -z_nd(i-1,j,k) - z_nd(i-1,j,k+1) );
243  return met_h_xi;
244 }

Referenced by ComputeStrain_T().

Here is the caller graph for this function:

◆ Compute_h_xi_AtIface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
112 {
113  amrex::Real met_h_xi = 0.125 * cellSizeInv[0] * (
114  (z_nd(i+1,j ,k ) - z_nd(i-1,j ,k )) + (z_nd(i+1,j ,k+1) - z_nd(i-1,j ,k+1)) +
115  (z_nd(i+1,j+1,k ) - z_nd(i-1,j+1,k )) + (z_nd(i+1,j+1,k+1) - z_nd(i-1,j+1,k+1)) );
116  return met_h_xi;
117 }

Referenced by DiffusionSrcForState_T(), erf_fast_rhs_MT(), erf_fast_rhs_T(), ERFPhysBCFunct_u::impose_vertical_xvel_bcs(), rotate_stress_tensor(), and MOSTAverage::set_rotated_fields().

Here is the caller graph for this function:

◆ Compute_h_xi_AtJface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtJface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
152 {
153  amrex::Real met_h_xi = 0.5 * cellSizeInv[0] * ( (z_nd(i+1,j,k ) - z_nd(i,j,k ))
154  + (z_nd(i+1,j,k+1) - z_nd(i,j,k+1)) );
155  return met_h_xi;
156 }

Referenced by ERFPhysBCFunct_v::impose_vertical_yvel_bcs().

Here is the caller graph for this function:

◆ Compute_h_xi_AtKface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtKface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
193 {
194  amrex::Real met_h_xi = 0.5 * cellSizeInv[0] * ( (z_nd(i+1,j ,k) - z_nd(i,j ,k))
195  + (z_nd(i+1,j+1,k) - z_nd(i,j+1,k)) );
196  return met_h_xi;
197 }

Referenced by DiffusionSrcForState_T().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtCellCenter()

AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real Compute_h_zeta_AtCellCenter ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
50 {
51  amrex::Real dzInv = cellSizeInv[2];
52  amrex::Real met_h_zeta = 0.25 * dzInv *
53  ( z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k+1) + z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
54  -z_nd(i+1,j,k ) - z_nd(i+1,j+1,k ) - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
55  return met_h_zeta;
56 }

Referenced by MOSTAverage::compute_region_averages(), ComputeDiffusivityMRF(), ComputeDiffusivityMYJ(), ComputeDiffusivityMYNN25(), ComputeDiffusivityMYNNEDMF(), ComputeDiffusivityYSU(), ComputeTurbulentViscosityLES(), ComputeTurbulentViscosityRANS(), and ERFPhysBCFunct_cons::impose_vertical_cons_bcs().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtEdgeCenterI()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterI ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
314 {
315  amrex::Real dzInv = cellSizeInv[2];
316  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1)
317  -z_nd(i,j,k-1) - z_nd(i+1,j,k-1) );
318  return met_h_zeta;
319 }

Referenced by AdvectionSrcForMom_TF(), AdvectionSrcForZMom(), ComputeStrain_T(), ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtEdgeCenterJ()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterJ ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
270 {
271  amrex::Real dzInv = cellSizeInv[2];
272  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
273  -z_nd(i,j,k-1) - z_nd(i,j+1,k-1) );
274  return met_h_zeta;
275 }

Referenced by AdvectionSrcForMom_TF(), AdvectionSrcForZMom(), ComputeStrain_T(), ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtEdgeCenterK()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterK ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
225 {
226  amrex::Real dzInv = cellSizeInv[2];
227  amrex::Real met_h_zeta = dzInv * (z_nd(i,j,k+1) - z_nd(i,j,k));
228  return met_h_zeta;
229 }

Referenced by AdvectionSrcForMom_TF(), AdvectionSrcForXMom(), AdvectionSrcForYMom(), ComputeStressConsVisc_T(), and ComputeStressVarVisc_T().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtIface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
99 {
100  amrex::Real met_h_zeta = 0.5 * cellSizeInv[2] * ( (z_nd(i,j ,k+1) - z_nd(i,j ,k))
101  + (z_nd(i,j+1,k+1) - z_nd(i,j+1,k)) );
102  return met_h_zeta;
103 }

Referenced by erf_fast_rhs_MT(), erf_fast_rhs_T(), erf_make_tau_terms(), erf_slow_rhs_pre(), ERFPhysBCFunct_u::impose_vertical_xvel_bcs(), and rotate_stress_tensor().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtJface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
139 {
140  amrex::Real met_h_zeta = 0.5 * cellSizeInv[2] * ( (z_nd(i ,j,k+1) - z_nd(i ,j,k ))
141  + (z_nd(i+1,j,k+1) - z_nd(i+1,j,k )) );
142  return met_h_zeta;
143 }

Referenced by erf_fast_rhs_MT(), erf_fast_rhs_T(), erf_make_tau_terms(), erf_slow_rhs_pre(), and ERFPhysBCFunct_v::impose_vertical_yvel_bcs().

Here is the caller graph for this function:

◆ Compute_h_zeta_AtKface()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface ( const int &  i,
const int &  j,
const int &  k,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  cellSizeInv,
const amrex::Array4< const amrex::Real > &  z_nd 
)
179 {
180  amrex::Real met_h_zeta = 0.125 * cellSizeInv[2] *
181  ( (z_nd(i ,j ,k+1) - z_nd(i ,j ,k-1)) + (z_nd(i+1,j ,k+1) - z_nd(i+1,j ,k-1)) +
182  (z_nd(i ,j+1,k+1) - z_nd( i,j+1,k-1)) + (z_nd(i+1,j+1,k+1) - z_nd(i+1,j+1,k-1)) );
183  return met_h_zeta;
184 }

Referenced by compute_gradp(), compute_gradp_interpz(), and DiffusionSrcForState_T().

Here is the caller graph for this function:

◆ Compute_Z_AtCellCenter()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtCellCenter ( const int &  i,
const int &  j,
const int &  k,
const amrex::Array4< const amrex::Real > &  z_nd 
)
355 {
356  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
357  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
358  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
359  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
360 
361  return z_cc;
362 }

Referenced by ERFPhysBCFunct_cons::impose_vertical_cons_bcs().

Here is the caller graph for this function:

◆ Compute_Z_AtWFace()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace ( const int &  i,
const int &  j,
const int &  k,
const amrex::Array4< const amrex::Real > &  z_nd 
)
370 {
371  const amrex::Real z_wf = 0.25*( z_nd(i ,j ,k ) + z_nd(i+1,j ,k )
372  + z_nd(i ,j+1,k ) + z_nd(i+1,j+1,k ) );
373 
374  return z_wf;
375 }

Referenced by DiffusionSrcForState_T(), OmegaFromW(), WFromOmega(), and ERF::Write2DPlotFile().

Here is the caller graph for this function:

◆ Compute_Zrel_AtCellCenter()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Zrel_AtCellCenter ( const int &  i,
const int &  j,
const int &  k,
const amrex::Array4< const amrex::Real > &  z_nd 
)
383 {
384  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
385  + z_nd(i+1,j ,k ) + z_nd(i+1,j ,k+1)
386  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
387  + z_nd(i+1,j+1,k ) + z_nd(i+1,j+1,k+1));
388 
389  // Note: we assume the z_nd array spans from the bottom to top of the domain
390  // i.e. no domain decomposition across processors in vertical direction
391  const amrex::Real z0_cc = 0.25*( z_nd(i ,j ,0) + z_nd(i ,j+1,0)
392  + z_nd(i+1,j ,0) + z_nd(i+1,j+1,0));
393 
394  return (z_cc - z0_cc);
395 }

Referenced by ComputeDiffusivityMRF(), ComputeDiffusivityMYJ(), ComputeDiffusivityMYNN25(), ComputeDiffusivityMYNNEDMF(), ComputeDiffusivityYSU(), SurfaceLayer::fill_qsurf_with_qsat(), ERFPhysBCFunct_cons::impose_vertical_cons_bcs(), and SurfaceLayer::init_tke_from_ustar().

Here is the caller graph for this function:

◆ init_default_zphys()

void init_default_zphys ( int  lev,
const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  z_phys_cc 
)

Routine to define default z_phys_nd and z_phys_cc

◆ init_which_terrain_grid()

void init_which_terrain_grid ( int  lev,
const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::Vector< amrex::Real > const &  z_levels_h 
)

◆ init_zlevels()

void init_zlevels ( amrex::Vector< amrex::Vector< amrex::Real >> &  zlevels_stag,
amrex::Vector< amrex::Vector< amrex::Real >> &  stretched_dz_h,
amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real >> &  stretched_dz_d,
amrex::Vector< amrex::Geometry > const &  geom,
amrex::Vector< amrex::IntVect > const &  ref_ratio,
const amrex::Real  grid_stretching_ratio,
const amrex::Real  zsurf,
const amrex::Real  dz0 
)

Utility routines for constructing terrain metric terms

◆ make_terrain_fitted_coords()

void make_terrain_fitted_coords ( int  lev,
const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::Vector< amrex::Real > const &  z_levels_h,
amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &  phys_bc_type 
)

◆ OmegaFromW()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real OmegaFromW ( int &  i,
int &  j,
int &  k,
amrex::Real  w,
const amrex::Array4< const amrex::Real > &  u_arr,
const amrex::Array4< const amrex::Real > &  v_arr,
const amrex::Array4< const amrex::Real > &  mf_u,
const amrex::Array4< const amrex::Real > &  mf_v,
const amrex::Array4< const amrex::Real > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv 
)
413 {
414  // This is dh/dxi at hi and lo edges
415  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
416  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
417  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
418  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
419  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
420  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
421 
422  // This is dh/deta at hi and lo edges
423  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
424  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
425  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
426  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
427  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
428  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
429 
430  // Use extrapolation instead of interpolation if at the bottom boundary
431  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
432  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
433  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
434  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
435  amrex::Real mf_u_hi = mf_u(i+1,j,0);
436  amrex::Real mf_u_lo = mf_u(i ,j,0);
437 
438  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
439  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
440  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
441  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
442  amrex::Real mf_v_hi = mf_v(i,j+1,0);
443  amrex::Real mf_v_lo = mf_v(i,j ,0);
444 
445  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
446  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
447 
448  amrex::Real omega = w - u_met - v_met;
449  return omega;
450 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:368
@ omega
Definition: ERF_Morrison.H:53

Referenced by erf_fast_rhs_MT(), erf_fast_rhs_T(), erf_make_tau_terms(), erf_slow_rhs_pre(), and ERF::project_momenta().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rotate_scalar_flux()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_scalar_flux ( const int &  i,
const int &  j,
const int &  klo,
const amrex::Real flux,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv,
const amrex::Array4< const amrex::Real > &  zphys_arr,
const amrex::Array4< amrex::Real > &  phi1_arr,
const amrex::Array4< amrex::Real > &  phi2_arr,
const amrex::Array4< amrex::Real > &  phi3_arr 
)
612 {
613  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
614  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
615  amrex::Real InvNorm = 1.0 / std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
616  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
617  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
618  phi3_arr(i,j,klo) = flux * InvNorm;
619 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:77
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:62

Referenced by SurfaceLayer::compute_SurfaceLayer_bcs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rotate_stress_tensor()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_stress_tensor ( const int &  i,
const int &  j,
const int &  klo,
const amrex::Real flux,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv,
const amrex::Array4< const amrex::Real > &  zphys_arr,
const amrex::Array4< const amrex::Real > &  u_arr,
const amrex::Array4< const amrex::Real > &  v_arr,
const amrex::Array4< const amrex::Real > &  w_arr,
const amrex::Array4< amrex::Real > &  tau11_arr,
const amrex::Array4< amrex::Real > &  tau22_arr,
const amrex::Array4< amrex::Real > &  tau33_arr,
const amrex::Array4< amrex::Real > &  tau12_arr,
const amrex::Array4< amrex::Real > &  tau21_arr,
const amrex::Array4< amrex::Real > &  tau13_arr,
const amrex::Array4< amrex::Real > &  tau31_arr,
const amrex::Array4< amrex::Real > &  tau23_arr,
const amrex::Array4< amrex::Real > &  tau32_arr 
)
642 {
643  // Unit-normal vector
644  amrex::Array1D<amrex::Real,0,2> n_hat;
645 
646  // Unit-tangent vectors
647  amrex::Array1D<amrex::Real,0,2> t_hat_1;
648  amrex::Array1D<amrex::Real,0,2> t_hat_2;
649  amrex::Array1D<amrex::Real,0,2> u_t_hat;
650 
651  // Final basis vector for right-hand coord sys
652  amrex::Array1D<amrex::Real,0,2> a_hat;
653 
654  // Rotation matrix
655  amrex::Array2D<amrex::Real,0,2,0,2> R_mat;
656 
657  // Metric data
658  amrex::Real h_xi = Compute_h_xi_AtIface(i, j, klo, dxInv, zphys_arr);
659  amrex::Real h_eta = Compute_h_eta_AtIface(i, j, klo, dxInv, zphys_arr);
660  amrex::Real h_zeta = Compute_h_zeta_AtIface(i, j, klo, dxInv, zphys_arr);
661 
662  // Populate the normal vector
663  amrex::Real Inormn = 1./std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
664  n_hat(0) = Inormn*h_xi; n_hat(1) = Inormn*h_eta; n_hat(2) = -Inormn;
665 
666  // Populate the tangent vectors
667  amrex::Real Inorm1 = 1./std::sqrt(1.0 + h_xi*h_xi);
668  amrex::Real Inorm2 = 1./std::sqrt(1.0 + h_eta*h_eta);
669  t_hat_1(0) = -Inorm1; t_hat_2(1) = -Inorm2;
670  t_hat_1(2) = -Inorm1*h_xi; t_hat_2(2) = -Inorm2*h_eta;
671 
672  // Populate the u_t vector
673  amrex::Real Norm_u_t = 0.0;
674  amrex::Real mag1 = (u_arr(i,j,klo) + h_xi *w_arr(i,j,klo))*Inorm1;
675  amrex::Real mag2 = (v_arr(i,j,klo) + h_eta*w_arr(i,j,klo))*Inorm2;
676  for (int icol(0); icol<3; ++icol) {
677  u_t_hat(icol) = mag1*t_hat_1(icol) + mag2*t_hat_2(icol);
678  Norm_u_t += u_t_hat(icol)*u_t_hat(icol);
679  }
680  for (int icol(0); icol<3; ++icol) {
681  u_t_hat(icol) /= std::sqrt(Norm_u_t);
682  }
683 
684  // Populate the a_hat vector
685  a_hat(0) = n_hat(1)*u_t_hat(2) - n_hat(2)*u_t_hat(1);
686  a_hat(1) = -(n_hat(0)*u_t_hat(2) - n_hat(2)*u_t_hat(0));
687  a_hat(2) = n_hat(0)*u_t_hat(1) - n_hat(1)*u_t_hat(0);
688 
689  // Copy column vectors into R_mat
690  int jrow;
691  jrow = 0;
692  for (int icol(0); icol<3; ++icol) {
693  R_mat(icol,jrow) = u_t_hat(icol);
694  }
695  jrow = 1;
696  for (int icol(0); icol<3; ++icol) {
697  R_mat(icol,jrow) = a_hat(icol);
698  }
699  jrow = 2;
700  for (int icol(0); icol<3; ++icol) {
701  R_mat(icol,jrow) = n_hat(icol);
702  }
703 
704  // Body-fixed tau
705  amrex::Real T11 = 2.0*R_mat(0,0)*R_mat(2,0)*flux;
706  amrex::Real T22 = 2.0*R_mat(0,1)*R_mat(2,1)*flux;
707  amrex::Real T33 = 2.0*R_mat(0,2)*R_mat(2,2)*flux;
708  amrex::Real T12 = (R_mat(0,0)*R_mat(2,1) + R_mat(2,0)*R_mat(0,1))*flux;
709  amrex::Real T13 = (R_mat(0,0)*R_mat(2,2) + R_mat(0,2)*R_mat(2,0))*flux;
710  amrex::Real T23 = (R_mat(0,1)*R_mat(2,2) + R_mat(0,2)*R_mat(2,1))*flux;
711 
712  // Rotated stress fluxes
713  tau11_arr(i,j,klo) = h_zeta*T11;
714  tau22_arr(i,j,klo) = h_zeta*T22;
715  tau33_arr(i,j,klo) = -h_xi*T13 - h_eta*T23 + T33;
716 
717  tau12_arr(i,j,klo) = h_zeta*T12;
718  tau21_arr(i,j,klo) = tau21_arr(i,j,klo);
719 
720  tau13_arr(i,j,klo) = -h_xi*T11 - h_eta*T12 + T13;
721  tau31_arr(i,j,klo) = h_zeta*T13;
722 
723  tau23_arr(i,j,klo) = -h_xi*T12 - h_eta*T22 + T23;
724  tau32_arr(i,j,klo) = h_zeta*T23;
725 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:109
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:96
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:123

Referenced by SurfaceLayer::compute_SurfaceLayer_bcs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ WFromOmega()

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega ( int &  i,
int &  j,
int &  k,
amrex::Real  omega,
const amrex::Array4< const amrex::Real > &  u_arr,
const amrex::Array4< const amrex::Real > &  v_arr,
const amrex::Array4< const amrex::Real > &  mf_u,
const amrex::Array4< const amrex::Real > &  mf_v,
const amrex::Array4< const amrex::Real > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv 
)
463 {
464  // This is dh/dxi at hi and lo edges
465  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
466  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
467  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
468  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
469  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
470  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
471 
472  // This is dh/deta at hi and lo edges
473  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
474  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
475  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
476  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
477  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
478  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
479 
480  // Use extrapolation instead of interpolation if at the bottom boundary
481  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
482  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
483  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
484  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
485  amrex::Real mf_u_hi = mf_u(i+1,j,0);
486  amrex::Real mf_u_lo = mf_u(i ,j,0);
487 
488  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
489  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
490  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
491  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
492  amrex::Real mf_v_hi = mf_v(i,j+1,0);
493  amrex::Real mf_v_lo = mf_v(i,j ,0);
494 
495  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
496  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
497 
498  amrex::Real w = omega + u_met + v_met;
499  return w;
500 }

Referenced by erf_fast_rhs_MT(), erf_fast_rhs_T(), ERFPhysBCFunct_w::impose_lateral_zvel_bcs(), ERFPhysBCFunct_w::impose_vertical_zvel_bcs(), and ERF::project_momenta().

Here is the call graph for this function:
Here is the caller graph for this function: