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_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 > 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, amrex::Real u, amrex::Real 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 > &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 
)
79 {
80  amrex::Real dyInv = cellSizeInv[1];
81  amrex::Real met_h_eta = 0.25 * dyInv *
82  ( 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)
83  -z_nd(i,j ,k) - z_nd(i+1,j ,k) - z_nd(i,j ,k+1) - z_nd(i+1,j ,k+1) );
84  return met_h_eta;
85 }

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

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

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

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

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 
)
64 {
65  amrex::Real dxInv = cellSizeInv[0];
66  amrex::Real met_h_xi = 0.25 * dxInv *
67  ( 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)
68  -z_nd(i ,j,k) - z_nd(i ,j+1,k) - z_nd(i ,j,k+1) - z_nd(i ,j+1,k+1) );
69  return met_h_xi;
70 }

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 
)
111 {
112  amrex::Real dxInv = cellSizeInv[0];
113  amrex::Real met_h_xi = 0.125 * dxInv *
114  ( 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)
115  -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) );
116  return met_h_xi;
117 }

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

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

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 
)
198 {
199  amrex::Real dxInv = cellSizeInv[0];
200  amrex::Real met_h_xi = 0.5 * dxInv * ( z_nd(i+1,j,k) + z_nd(i+1,j+1,k)
201  -z_nd(i ,j,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 
)
49 {
50  amrex::Real dzInv = cellSizeInv[2];
51  amrex::Real met_h_zeta = 0.25 * dzInv *
52  ( 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)
53  -z_nd(i+1,j,k ) - z_nd(i+1,j+1,k ) - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
54  return met_h_zeta;
55 }

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

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

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

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

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 
)
183 {
184  amrex::Real dzInv = cellSizeInv[2];
185  amrex::Real met_h_zeta = 0.125 * dzInv *
186  ( 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)
187  -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) );
188  return met_h_zeta;
189 }

Referenced by DiffusionSrcForState_T(), and erf_slow_rhs_pre().

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 
)
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  // Note: we assume the z_nd array spans from the bottom to top of the domain
368  // i.e. no domain decomposition across processors in vertical direction
369  const amrex::Real z0_cc = 0.25*( z_nd(i ,j ,0) + z_nd(i ,j+1,0)
370  + z_nd(i+1,j ,0) + z_nd(i+1,j+1,0));
371 
372  return (z_cc - z0_cc);
373 }

Referenced by ComputeDiffusivityMYNN25(), ComputeDiffusivityMYNNEDMF(), and ComputeDiffusivityYSU().

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 >  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv 
)

Define omega given u,v and w

386 {
387  // This is dh/dxi at z-face (i,j,k-1/2)
388  amrex::Real met_zlo_xi = 0.5 * dxInv[0] *
389  ( z_nd(i+1,j+1,k ) + z_nd(i+1,j ,k ) // hi i, lo k
390  -z_nd(i ,j+1,k ) - z_nd(i ,j ,k ) ); // lo i, lo k
391 
392  // This is dh/deta at z-face (i,j,k-1/2)
393  amrex::Real met_zlo_eta = 0.5 * dxInv[1] *
394  ( z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k ) // hi j, lo k
395  -z_nd(i+1,j ,k ) - z_nd(i ,j ,k ) ); // lo j, lo k
396 
397  // Slip BC or moving terrain
398  // Use extrapolation instead of interpolation if at the bottom boundary
399  amrex::Real u = (k == 0) ? 1.5 * (0.5*(u_arr(i,j,k)+u_arr(i+1,j,k))) - 0.5*(0.5*(u_arr(i,j,k+1)+u_arr(i+1,j,k+1))) :
400  0.25 * ( u_arr(i,j,k-1) + u_arr(i+1,j,k-1) + u_arr(i,j,k) + u_arr(i+1,j,k) );
401  amrex::Real v = (k == 0) ? 1.5 * (0.5*(v_arr(i,j,k)+v_arr(i,j+1,k))) - 0.5*(0.5*(v_arr(i,j,k+1)+v_arr(i,j+1,k+1))) :
402  0.25 * ( v_arr(i,j,k-1) + v_arr(i,j+1,k-1) + v_arr(i,j,k) + v_arr(i,j+1,k) );
403 
404  amrex::Real omega = w - met_zlo_xi * u - met_zlo_eta * v;
405  return omega;
406 }
@ omega
Definition: ERF_Morrison.H:45

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

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 
)
470 {
471  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
472  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
473  amrex::Real InvNorm = 1.0 / std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
474  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
475  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
476  phi3_arr(i,j,klo) = flux * InvNorm;
477 }
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:76
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:61

Referenced by ABLMost::compute_most_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 
)
500 {
501  // Unit-normal vector
502  amrex::Array1D<amrex::Real,0,2> n_hat;
503 
504  // Unit-tangent vectors
505  amrex::Array1D<amrex::Real,0,2> t_hat_1;
506  amrex::Array1D<amrex::Real,0,2> t_hat_2;
507  amrex::Array1D<amrex::Real,0,2> u_t_hat;
508 
509  // Final basis vector for right-hand coord sys
510  amrex::Array1D<amrex::Real,0,2> a_hat;
511 
512  // Rotation matrix
513  amrex::Array2D<amrex::Real,0,2,0,2> R_mat;
514 
515  // Metric data
516  amrex::Real h_xi = Compute_h_xi_AtIface(i, j, klo, dxInv, zphys_arr);
517  amrex::Real h_eta = Compute_h_eta_AtIface(i, j, klo, dxInv, zphys_arr);
518 
519  // Populate the normal vector
520  amrex::Real Inormn = 1./std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
521  n_hat(0) = -Inormn*h_xi; n_hat(1) = -Inormn*h_eta; n_hat(2) = Inormn;
522 
523  // Populate the tangent vectors
524  amrex::Real Inorm1 = 1./std::sqrt(1.0 + h_xi*h_xi);
525  amrex::Real Inorm2 = 1./std::sqrt(1.0 + h_eta*h_eta);
526  t_hat_1(0) = Inorm1; t_hat_2(1) = Inorm2;
527  t_hat_1(2) = Inorm1*h_xi; t_hat_2(2) = Inorm2*h_eta;
528 
529  // Populate the u_t vector
530  amrex::Real Norm_u_t = 0.0;
531  amrex::Real mag1 = (u_arr(i,j,klo) + h_xi *w_arr(i,j,klo))*Inorm1;
532  amrex::Real mag2 = (v_arr(i,j,klo) + h_eta*w_arr(i,j,klo))*Inorm2;
533  for (int icol(0); icol<3; ++icol) {
534  u_t_hat(icol) = mag1*t_hat_1(icol) + mag2*t_hat_2(icol);
535  Norm_u_t += u_t_hat(icol)*u_t_hat(icol);
536  }
537  for (int icol(0); icol<3; ++icol) {
538  u_t_hat(icol) /= std::sqrt(Norm_u_t);
539  }
540 
541  // Populate the a_hat vector
542  a_hat(0) = n_hat(1)*u_t_hat(2) - n_hat(2)*u_t_hat(1);
543  a_hat(1) = -(n_hat(0)*u_t_hat(2) - n_hat(2)*u_t_hat(0));
544  a_hat(2) = n_hat(1)*u_t_hat(1) - n_hat(1)*u_t_hat(0);
545 
546  // Copy column vectors into R_mat
547  int jrow;
548  jrow = 0;
549  for (int icol(0); icol<3; ++icol) {
550  R_mat(icol,jrow) = u_t_hat(icol);
551  }
552  jrow = 1;
553  for (int icol(0); icol<3; ++icol) {
554  R_mat(icol,jrow) = a_hat(icol);
555  }
556  jrow = 2;
557  for (int icol(0); icol<3; ++icol) {
558  R_mat(icol,jrow) = n_hat(icol);
559  }
560 
561  // Assign the stresses
562  tau11_arr(i,j,klo) = 2.0*R_mat(0,0)*R_mat(2,0)*flux;
563  tau22_arr(i,j,klo) = 2.0*R_mat(0,1)*R_mat(2,1)*flux;
564  tau33_arr(i,j,klo) = 2.0*R_mat(0,2)*R_mat(2,2)*flux;
565 
566  tau12_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,1) + R_mat(2,0)*R_mat(0,1))*flux;
567  tau21_arr(i,j,klo) = tau12_arr(i,j,klo);
568 
569  tau13_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,2) + R_mat(0,2)*R_mat(2,0))*flux;
570  tau31_arr(i,j,klo) = tau13_arr(i,j,klo);
571 
572  tau23_arr(i,j,klo) = (R_mat(0,1)*R_mat(2,2) + R_mat(0,2)*R_mat(2,1))*flux;
573  tau32_arr(i,j,klo) = tau23_arr(i,j,klo);
574 }
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:108
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 ABLMost::compute_most_bcs().

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

◆ WFromOmega() [1/2]

AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega ( int  i,
int  j,
int  k,
amrex::Real  omega,
amrex::Real  u,
amrex::Real  v,
const amrex::Array4< const amrex::Real > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv 
)

Define w given scalar u,v and omega

418 {
419  // This is dh/dxi at z-face (i,j,k-1/2)
420  amrex::Real met_zlo_xi = 0.5 * dxInv[0] *
421  ( z_nd(i+1,j+1,k ) + z_nd(i+1,j ,k ) // hi i, lo k
422  -z_nd(i ,j+1,k ) - z_nd(i ,j ,k ) ); // lo i, lo k
423 
424  // This is dh/deta at z-face (i,j,k-1/2)
425  amrex::Real met_zlo_eta = 0.5 * dxInv[1] *
426  ( z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k ) // hi j, lo k
427  -z_nd(i+1,j ,k ) - z_nd(i ,j ,k ) ); // lo j, lo k
428 
429  amrex::Real w = omega + met_zlo_xi * u + met_zlo_eta * v;
430  return w;
431 }

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

Here is the caller graph for this function:

◆ WFromOmega() [2/2]

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 > &  z_nd,
const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &  dxInv 
)

Define w given u and v arrays and scalar omega

444 {
445  // Use extrapolation instead of interpolation if at the bottom boundary
446  amrex::Real u = (k == 0) ? 1.5 * (0.5*(u_arr(i,j,k)+u_arr(i+1,j,k))) - 0.5*(0.5*(u_arr(i,j,k+1)+u_arr(i+1,j,k+1))) :
447  0.25 * ( u_arr(i,j,k-1) + u_arr(i+1,j,k-1) + u_arr(i,j,k) + u_arr(i+1,j,k) );
448  amrex::Real v = (k == 0) ? 1.5 * (0.5*(v_arr(i,j,k)+v_arr(i,j+1,k))) - 0.5*(0.5*(v_arr(i,j,k+1)+v_arr(i,j+1,k+1))) :
449  0.25 * ( v_arr(i,j,k-1) + v_arr(i,j+1,k-1) + v_arr(i,j,k) + v_arr(i,j+1,k) );
450 
451  amrex::Real w = WFromOmega(i,j,k,omega,u,v,z_nd,dxInv);
452  return w;
453 }
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int i, int j, int k, amrex::Real omega, amrex::Real u, amrex::Real v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:414
Here is the call graph for this function: