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::Real get_dzmin_terrain (amrex::MultiFab &z_phys_nd)
 
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 
)
86 {
87  amrex::Real dyInv = cellSizeInv[1];
88  amrex::Real met_h_eta = 0.25 * dyInv *
89  ( 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)
90  -z_nd(i,j ,k) - z_nd(i+1,j ,k) - z_nd(i,j ,k+1) - z_nd(i+1,j ,k+1) );
91  return met_h_eta;
92 }
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 
)
347 {
348  amrex::Real dyInv = cellSizeInv[1];
349  amrex::Real met_h_eta = 0.25 * dyInv *
350  ( z_nd(i+1,j+1,k) + z_nd(i,j+1,k)
351  -z_nd(i+1,j-1,k) - z_nd(i,j-1,k) );
352  return met_h_eta;
353 }

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

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

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

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

Referenced by DiffusionSrcForState_T(), erf_substep_MT(), erf_substep_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 
)
212 {
213  amrex::Real met_h_eta = 0.5 * cellSizeInv[1] * ( (z_nd(i ,j+1,k) - z_nd(i ,j,k))
214  + (z_nd(i+1,j+1,k) - z_nd(i+1,j,k)) );
215  return met_h_eta;
216 }

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 
)
71 {
72  amrex::Real dxInv = cellSizeInv[0];
73  amrex::Real met_h_xi = 0.25 * dxInv *
74  ( 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)
75  -z_nd(i ,j,k) - z_nd(i ,j+1,k) - z_nd(i ,j,k+1) - z_nd(i ,j+1,k+1) );
76  return met_h_xi;
77 }

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

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

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

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 
)
118 {
119  amrex::Real met_h_xi = 0.125 * cellSizeInv[0] * (
120  (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)) +
121  (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)) );
122  return met_h_xi;
123 }

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

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

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 
)
56 {
57  amrex::Real dzInv = cellSizeInv[2];
58  amrex::Real met_h_zeta = 0.25 * dzInv *
59  ( 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)
60  -z_nd(i+1,j,k ) - z_nd(i+1,j+1,k ) - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
61  return met_h_zeta;
62 }

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

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

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

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 
)
105 {
106  amrex::Real met_h_zeta = 0.5 * cellSizeInv[2] * ( (z_nd(i,j ,k+1) - z_nd(i,j ,k))
107  + (z_nd(i,j+1,k+1) - z_nd(i,j+1,k)) );
108  return met_h_zeta;
109 }

Referenced by erf_make_tau_terms(), erf_slow_rhs_pre(), erf_substep_MT(), erf_substep_T(), 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 
)
145 {
146  amrex::Real met_h_zeta = 0.5 * cellSizeInv[2] * ( (z_nd(i ,j,k+1) - z_nd(i ,j,k ))
147  + (z_nd(i+1,j,k+1) - z_nd(i+1,j,k )) );
148  return met_h_zeta;
149 }

Referenced by erf_make_tau_terms(), erf_slow_rhs_pre(), erf_substep_MT(), erf_substep_T(), 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 met_h_zeta = 0.125 * cellSizeInv[2] *
187  ( (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)) +
188  (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)) );
189  return met_h_zeta;
190 }

Referenced by compute_gradp(), compute_gradp_interpz(), DiffusionSrcForState_T(), and ImplicitDiffForState_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 
)
361 {
362  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
363  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
364  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
365  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
366 
367  return z_cc;
368 }

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

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

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

Here is the caller graph for this function:

◆ get_dzmin_terrain()

amrex::Real get_dzmin_terrain ( amrex::MultiFab &  z_phys_nd)

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

Referenced by erf_make_tau_terms(), erf_slow_rhs_pre(), erf_substep_MT(), erf_substep_T(), 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 
)
618 {
619  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
620  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
621  amrex::Real InvNorm = 1.0 / std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
622  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
623  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
624  phi3_arr(i,j,klo) = flux * InvNorm;
625 }
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:83
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:68

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

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

Referenced by erf_substep_MT(), erf_substep_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: