ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 }

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 
)
350 {
351  amrex::Real dyInv = cellSizeInv[1];
352  amrex::Real met_h_eta = 0.25 * dyInv *
353  ( z_nd(i+1,j+1,k) + z_nd(i,j+1,k)
354  -z_nd(i+1,j-1,k) - z_nd(i,j-1,k) );
355  return met_h_eta;
356 }

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 
)
308 {
309  amrex::Real dyInv = cellSizeInv[1];
310  amrex::Real met_h_eta = dyInv * ( z_nd(i,j+1,k) - z_nd(i,j,k) );
311  return met_h_eta;
312 }

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 
)
262 {
263  amrex::Real dyInv = cellSizeInv[1];
264  amrex::Real met_h_eta = 0.25 * dyInv *
265  ( z_nd(i,j+1,k) + z_nd(i,j+1,k+1)
266  -z_nd(i,j-1,k) - z_nd(i,j-1,k+1) );
267  return met_h_eta;
268 }

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 
)
128 {
129  amrex::Real dyInv = cellSizeInv[1];
130  amrex::Real met_h_eta = 0.5 * dyInv * ( z_nd(i,j+1,k ) + z_nd(i,j+1,k+1)
131  - z_nd(i,j ,k ) - z_nd(i,j ,k+1) );
132  return met_h_eta;
133 }

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 
)
170 {
171  amrex::Real dyInv = cellSizeInv[1];
172  amrex::Real met_h_eta = 0.125 * dyInv *
173  ( z_nd(i,j+1,k) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k) + z_nd(i+1,j+1,k+1)
174  -z_nd(i,j-1,k) - z_nd(i,j-1,k+1) - z_nd(i+1,j-1,k) - z_nd(i+1,j-1,k+1) );
175  return met_h_eta;
176 }

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 
)
214 {
215  amrex::Real dyInv = cellSizeInv[1];
216  amrex::Real met_h_eta = 0.5 * dyInv * ( z_nd(i,j+1,k) + z_nd(i+1,j+1,k)
217  -z_nd(i,j ,k) - z_nd(i+1,j ,k) );
218  return met_h_eta;
219 }

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 
)
337 {
338  amrex::Real dxInv = cellSizeInv[0];
339  amrex::Real met_h_xi = dxInv * ( z_nd(i+1,j,k) - z_nd(i,j,k) );
340  return met_h_xi;
341 }

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 
)
293 {
294  amrex::Real dxInv = cellSizeInv[0];
295  amrex::Real met_h_xi = 0.25 * dxInv *
296  ( z_nd(i+1,j+1,k) + z_nd(i+1,j ,k)
297  -z_nd(i-1,j+1,k) - z_nd(i-1,j ,k) );
298  return met_h_xi;
299 }

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 
)
247 {
248  amrex::Real dxInv = cellSizeInv[0];
249  amrex::Real met_h_xi = 0.25 * dxInv *
250  ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1)
251  -z_nd(i-1,j,k) - z_nd(i-1,j,k+1) );
252  return met_h_xi;
253 }

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 
)
113 {
114  amrex::Real dxInv = cellSizeInv[0];
115  amrex::Real met_h_xi = 0.125 * dxInv *
116  ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k) + z_nd(i+1,j+1,k+1)
117  -z_nd(i-1,j,k) - z_nd(i-1,j,k+1) - z_nd(i-1,j+1,k) - z_nd(i-1,j+1,k+1) );
118  return met_h_xi;
119 }

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 
)
156 {
157  amrex::Real dxInv = cellSizeInv[0];
158  amrex::Real met_h_xi = 0.5 * dxInv * ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1)
159  -z_nd(i ,j,k) - z_nd(i ,j,k+1) );
160  return met_h_xi;
161 }

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 
)
200 {
201  amrex::Real dxInv = cellSizeInv[0];
202  amrex::Real met_h_xi = 0.5 * dxInv * ( z_nd(i+1,j,k) + z_nd(i+1,j+1,k)
203  -z_nd(i ,j,k) - z_nd(i ,j+1,k) );
204  return met_h_xi;
205 }

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(), 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 
)
323 {
324  amrex::Real dzInv = cellSizeInv[2];
325  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1)
326  -z_nd(i,j,k-1) - z_nd(i+1,j,k-1) );
327  return met_h_zeta;
328 }

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 
)
279 {
280  amrex::Real dzInv = cellSizeInv[2];
281  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
282  -z_nd(i,j,k-1) - z_nd(i,j+1,k-1) );
283  return met_h_zeta;
284 }

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 
)
234 {
235  amrex::Real dzInv = cellSizeInv[2];
236  amrex::Real met_h_zeta = dzInv * (z_nd(i,j,k+1) - z_nd(i,j,k));
237  return met_h_zeta;
238 }

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 dzInv = cellSizeInv[2];
101  amrex::Real met_h_zeta = 0.5 * dzInv * ( z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
102  - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
103  return met_h_zeta;
104 }

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

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 
)
142 {
143  amrex::Real dzInv = cellSizeInv[2];
144  amrex::Real met_h_zeta = 0.5 * dzInv * ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1)
145  - z_nd(i,j,k ) - z_nd(i+1,j,k ) );
146  return met_h_zeta;
147 }

Referenced by DiffusionSrcForState_T(), 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 
)
185 {
186  amrex::Real dzInv = cellSizeInv[2];
187  amrex::Real met_h_zeta = 0.125 * dzInv *
188  ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1)
189  -z_nd(i,j,k-1) - z_nd(i+1,j,k-1) - z_nd(i,j+1,k-1) - z_nd(i+1,j+1,k-1) );
190  return met_h_zeta;
191 }

Referenced by compute_gradp(), 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 
)
364 {
365  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
366  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
367  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
368  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
369 
370  return z_cc;
371 }

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 
)
379 {
380  const amrex::Real z_wf = 0.25*( z_nd(i ,j ,k ) + z_nd(i+1,j ,k )
381  + z_nd(i ,j+1,k ) + z_nd(i+1,j+1,k ) );
382 
383  return z_wf;
384 }

Referenced by OmegaFromW(), and WFromOmega().

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 
)
392 {
393  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
394  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
395  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
396  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
397 
398  // Note: we assume the z_nd array spans from the bottom to top of the domain
399  // i.e. no domain decomposition across processors in vertical direction
400  const amrex::Real z0_cc = 0.25*( z_nd(i ,j ,0) + z_nd(i ,j+1,0)
401  + z_nd(i+1,j ,0) + z_nd(i+1,j+1,0));
402 
403  return (z_cc - z0_cc);
404 }

Referenced by ComputeDiffusivityMYNN25(), ComputeDiffusivityMYNNEDMF(), ComputeDiffusivityYSU(), and ERFPhysBCFunct_cons::impose_vertical_cons_bcs().

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 
)
422 {
423  // This is dh/dxi at hi and lo edges
424  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
425  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
426  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
427  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
428  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
429  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
430 
431  // This is dh/deta at hi and lo edges
432  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
433  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
434  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
435  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
436  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
437  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
438 
439  // Use extrapolation instead of interpolation if at the bottom boundary
440  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
441  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
442  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
443  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
444  amrex::Real mf_u_hi = mf_u(i+1,j,0);
445  amrex::Real mf_u_lo = mf_u(i ,j,0);
446 
447  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
448  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
449  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
450  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
451  amrex::Real mf_v_hi = mf_v(i,j+1,0);
452  amrex::Real mf_v_lo = mf_v(i,j ,0);
453 
454  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
455  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
456 
457  amrex::Real omega = w - u_met - v_met;
458  return omega;
459 }
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:377
@ omega
Definition: ERF_Morrison.H:53

Referenced by erf_fast_rhs_MT(), erf_fast_rhs_T(), erf_make_tau_terms(), and erf_slow_rhs_pre().

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 
)
621 {
622  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
623  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
624  amrex::Real InvNorm = 1.0 / std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
625  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
626  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
627  phi3_arr(i,j,klo) = flux * InvNorm;
628 }
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 
)
651 {
652  // Unit-normal vector
653  amrex::Array1D<amrex::Real,0,2> n_hat;
654 
655  // Unit-tangent vectors
656  amrex::Array1D<amrex::Real,0,2> t_hat_1;
657  amrex::Array1D<amrex::Real,0,2> t_hat_2;
658  amrex::Array1D<amrex::Real,0,2> u_t_hat;
659 
660  // Final basis vector for right-hand coord sys
661  amrex::Array1D<amrex::Real,0,2> a_hat;
662 
663  // Rotation matrix
664  amrex::Array2D<amrex::Real,0,2,0,2> R_mat;
665 
666  // Metric data
667  amrex::Real h_xi = Compute_h_xi_AtIface(i, j, klo, dxInv, zphys_arr);
668  amrex::Real h_eta = Compute_h_eta_AtIface(i, j, klo, dxInv, zphys_arr);
669 
670  // Populate the normal vector
671  amrex::Real Inormn = 1./std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
672  n_hat(0) = -Inormn*h_xi; n_hat(1) = -Inormn*h_eta; n_hat(2) = Inormn;
673 
674  // Populate the tangent vectors
675  amrex::Real Inorm1 = 1./std::sqrt(1.0 + h_xi*h_xi);
676  amrex::Real Inorm2 = 1./std::sqrt(1.0 + h_eta*h_eta);
677  t_hat_1(0) = Inorm1; t_hat_2(1) = Inorm2;
678  t_hat_1(2) = Inorm1*h_xi; t_hat_2(2) = Inorm2*h_eta;
679 
680  // Populate the u_t vector
681  amrex::Real Norm_u_t = 0.0;
682  amrex::Real mag1 = (u_arr(i,j,klo) + h_xi *w_arr(i,j,klo))*Inorm1;
683  amrex::Real mag2 = (v_arr(i,j,klo) + h_eta*w_arr(i,j,klo))*Inorm2;
684  for (int icol(0); icol<3; ++icol) {
685  u_t_hat(icol) = mag1*t_hat_1(icol) + mag2*t_hat_2(icol);
686  Norm_u_t += u_t_hat(icol)*u_t_hat(icol);
687  }
688  for (int icol(0); icol<3; ++icol) {
689  u_t_hat(icol) /= std::sqrt(Norm_u_t);
690  }
691 
692  // Populate the a_hat vector
693  a_hat(0) = n_hat(1)*u_t_hat(2) - n_hat(2)*u_t_hat(1);
694  a_hat(1) = -(n_hat(0)*u_t_hat(2) - n_hat(2)*u_t_hat(0));
695  a_hat(2) = n_hat(1)*u_t_hat(1) - n_hat(1)*u_t_hat(0);
696 
697  // Copy column vectors into R_mat
698  int jrow;
699  jrow = 0;
700  for (int icol(0); icol<3; ++icol) {
701  R_mat(icol,jrow) = u_t_hat(icol);
702  }
703  jrow = 1;
704  for (int icol(0); icol<3; ++icol) {
705  R_mat(icol,jrow) = a_hat(icol);
706  }
707  jrow = 2;
708  for (int icol(0); icol<3; ++icol) {
709  R_mat(icol,jrow) = n_hat(icol);
710  }
711 
712  // Assign the stresses
713  tau11_arr(i,j,klo) = 2.0*R_mat(0,0)*R_mat(2,0)*flux;
714  tau22_arr(i,j,klo) = 2.0*R_mat(0,1)*R_mat(2,1)*flux;
715  tau33_arr(i,j,klo) = 2.0*R_mat(0,2)*R_mat(2,2)*flux;
716 
717  tau12_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,1) + R_mat(2,0)*R_mat(0,1))*flux;
718  tau21_arr(i,j,klo) = tau12_arr(i,j,klo);
719 
720  tau13_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,2) + R_mat(0,2)*R_mat(2,0))*flux;
721  tau31_arr(i,j,klo) = tau13_arr(i,j,klo);
722 
723  tau23_arr(i,j,klo) = (R_mat(0,1)*R_mat(2,2) + R_mat(0,2)*R_mat(2,1))*flux;
724  tau32_arr(i,j,klo) = tau23_arr(i,j,klo);
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:110
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:125

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 
)
472 {
473  // This is dh/dxi at hi and lo edges
474  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
475  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
476  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
477  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
478  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
479  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
480 
481  // This is dh/deta at hi and lo edges
482  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
483  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
484  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
485  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
486  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
487  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
488 
489  // Use extrapolation instead of interpolation if at the bottom boundary
490  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
491  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
492  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
493  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
494  amrex::Real mf_u_hi = mf_u(i+1,j,0);
495  amrex::Real mf_u_lo = mf_u(i ,j,0);
496 
497  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
498  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
499  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
500  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
501  amrex::Real mf_v_hi = mf_v(i,j+1,0);
502  amrex::Real mf_v_lo = mf_v(i,j ,0);
503 
504  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
505  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
506 
507  amrex::Real w = omega + u_met + v_met;
508  return w;
509 }

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_velocities().

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