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 <ERF_Constants.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, amrex::Real z_offset=zero)
 
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 
)
88 {
89  amrex::Real dyInv = cellSizeInv[1];
90  amrex::Real met_h_eta = fourth * dyInv *
91  ( 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)
92  -z_nd(i,j ,k) - z_nd(i+1,j ,k) - z_nd(i,j ,k+1) - z_nd(i+1,j ,k+1) );
93  return met_h_eta;
94 }
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
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 
)
349 {
350  amrex::Real dyInv = cellSizeInv[1];
351  amrex::Real met_h_eta = fourth * dyInv *
352  ( z_nd(i+1,j+1,k) + z_nd(i,j+1,k)
353  -z_nd(i+1,j-1,k) - z_nd(i,j-1,k) );
354  return met_h_eta;
355 }

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

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

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 
)
134 {
135  amrex::Real met_h_eta = myhalf * cellSizeInv[1] * ( (z_nd(i,j+1,k ) - z_nd(i,j,k ))
136  + (z_nd(i,j+1,k+1) - z_nd(i,j,k+1)) );
137  return met_h_eta;
138 }
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11

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

Referenced by DiffusionSrcForState_T(), erf_substep_MT(), erf_substep_T(), ERFPhysBCFunct_v::impose_vertical_yvel_bcs(), ERF::poisson_wall_dist(), 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 met_h_eta = myhalf * cellSizeInv[1] * ( (z_nd(i ,j+1,k) - z_nd(i ,j,k))
216  + (z_nd(i+1,j+1,k) - z_nd(i+1,j,k)) );
217  return met_h_eta;
218 }

Referenced by DiffusionSrcForState_T(), and ERF::poisson_wall_dist().

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 
)
73 {
74  amrex::Real dxInv = cellSizeInv[0];
75  amrex::Real met_h_xi = fourth * dxInv *
76  ( 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)
77  -z_nd(i ,j,k) - z_nd(i ,j+1,k) - z_nd(i ,j,k+1) - z_nd(i ,j+1,k+1) );
78  return met_h_xi;
79 }
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17

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

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

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

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

Referenced by DiffusionSrcForState_T(), erf_substep_MT(), erf_substep_T(), ERFPhysBCFunct_u::impose_vertical_xvel_bcs(), ERF::poisson_wall_dist(), 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 
)
160 {
161  amrex::Real met_h_xi = myhalf * cellSizeInv[0] * ( (z_nd(i+1,j,k ) - z_nd(i,j,k ))
162  + (z_nd(i+1,j,k+1) - z_nd(i,j,k+1)) );
163  return met_h_xi;
164 }

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

Referenced by DiffusionSrcForState_T(), and ERF::poisson_wall_dist().

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

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

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

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

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

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

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

Referenced by erf_make_tau_terms(), erf_slow_rhs_pre(), erf_substep_MT(), erf_substep_T(), ERFPhysBCFunct_v::impose_vertical_yvel_bcs(), and ERF::poisson_wall_dist().

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 
)
187 {
188  amrex::Real met_h_zeta = amrex::Real(0.125) * cellSizeInv[2] *
189  ( (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)) +
190  (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)) );
191  return met_h_zeta;
192 }

Referenced by compute_gradp(), compute_gradp_interpz(), DiffusionSrcForState_T(), ImplicitDiffForMomLU_T(), ImplicitDiffForStateLU_T(), and ERF::poisson_wall_dist().

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

Referenced by ERFPhysBCFunct_cons::impose_vertical_cons_bcs(), and SurfaceLayer::init_tke_from_ustar().

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

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

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

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,
amrex::Real  z_offset = zero 
)

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 
)
421 {
422  // This is dh/dxi at hi and lo edges
423  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
424  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
425  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
426  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
427  amrex::Real met_xi = (amrex::Real(1.0/amrex::Real(12.0))*(z_x_m2 - z_x_p2)
428  + amrex::Real(8.0/amrex::Real(12.0))*(z_x_p1 - z_x_m1)) * dxInv[0];
429 
430  // This is dh/deta at hi and lo edges
431  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
432  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
433  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
434  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
435  amrex::Real met_eta = (amrex::Real(1.0/amrex::Real(12.0))*(z_y_m2 - z_y_p2)
436  + amrex::Real(8.0/amrex::Real(12.0))*(z_y_p1 - z_y_m1)) * dxInv[1];
437 
438  // Use extrapolation instead of interpolation if at the bottom boundary
439  amrex::Real u_hi = (k == 0) ? amrex::Real(1.5) * u_arr(i+1,j ,k ) - myhalf * u_arr(i+1,j ,k+1) :
440  myhalf * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
441  amrex::Real u_lo = (k == 0) ? amrex::Real(1.5) * u_arr(i ,j ,k ) - myhalf * u_arr(i ,j ,k+1) :
442  myhalf * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
443  amrex::Real mf_u_hi = mf_u(i+1,j,0);
444  amrex::Real mf_u_lo = mf_u(i ,j,0);
445 
446  amrex::Real v_hi = (k == 0) ? amrex::Real(1.5) * v_arr(i ,j+1,k ) - myhalf * v_arr(i ,j+1,k+1) :
447  myhalf * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
448  amrex::Real v_lo = (k == 0) ? amrex::Real(1.5) * v_arr(i ,j ,k ) - myhalf * v_arr(i ,j ,k+1) :
449  myhalf * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
450  amrex::Real mf_v_hi = mf_v(i,j+1,0);
451  amrex::Real mf_v_lo = mf_v(i,j ,0);
452 
453  amrex::Real u_met = met_xi * myhalf * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
454  amrex::Real v_met = met_eta * myhalf * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
455 
456  amrex::Real omega = w - u_met - v_met;
457  return omega;
458 }
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:376
@ 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 
)
620 {
621  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
622  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
623  amrex::Real InvNorm = one / std::sqrt(one + h_xi*h_xi + h_eta*h_eta);
624  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
625  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
626  phi3_arr(i,j,klo) = flux * InvNorm;
627 }
constexpr amrex::Real one
Definition: ERF_Constants.H:7
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:85
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:70

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 
)
650 {
651  // Unit-normal vector
652  amrex::Array1D<amrex::Real,0,2> n_hat;
653 
654  // Unit-tangent vectors
655  amrex::Array1D<amrex::Real,0,2> t_hat_1;
656  amrex::Array1D<amrex::Real,0,2> t_hat_2;
657  amrex::Array1D<amrex::Real,0,2> u_t_hat;
658 
659  // Final basis vector for right-hand coord sys
660  amrex::Array1D<amrex::Real,0,2> a_hat;
661 
662  // Rotation matrix
663  amrex::Array2D<amrex::Real,0,2,0,2> R_mat;
664 
665  // Metric data
666  amrex::Real h_xi = Compute_h_xi_AtIface(i, j, klo, dxInv, zphys_arr);
667  amrex::Real h_eta = Compute_h_eta_AtIface(i, j, klo, dxInv, zphys_arr);
668  amrex::Real h_zeta = Compute_h_zeta_AtIface(i, j, klo, dxInv, zphys_arr);
669 
670  // Populate the normal vector
671  amrex::Real Inormn = one/std::sqrt(one + 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 = one/std::sqrt(one + h_xi*h_xi);
676  amrex::Real Inorm2 = one/std::sqrt(one + 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 = zero;
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(0)*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  // Body-fixed tau
713  amrex::Real T11 = two*R_mat(0,0)*R_mat(2,0)*flux;
714  amrex::Real T22 = two*R_mat(0,1)*R_mat(2,1)*flux;
715  amrex::Real T33 = two*R_mat(0,2)*R_mat(2,2)*flux;
716  amrex::Real T12 = (R_mat(0,0)*R_mat(2,1) + R_mat(2,0)*R_mat(0,1))*flux;
717  amrex::Real T13 = (R_mat(0,0)*R_mat(2,2) + R_mat(0,2)*R_mat(2,0))*flux;
718  amrex::Real T23 = (R_mat(0,1)*R_mat(2,2) + R_mat(0,2)*R_mat(2,1))*flux;
719 
720  // Rotated stress fluxes
721  tau11_arr(i,j,klo) = h_zeta*T11;
722  tau22_arr(i,j,klo) = h_zeta*T22;
723  tau33_arr(i,j,klo) = -h_xi*T13 - h_eta*T23 + T33;
724 
725  tau12_arr(i,j,klo) = h_zeta*T12;
726  tau21_arr(i,j,klo) = tau21_arr(i,j,klo);
727 
728  tau13_arr(i,j,klo) = -h_xi*T11 - h_eta*T12 + T13;
729  tau31_arr(i,j,klo) = h_zeta*T13;
730 
731  tau23_arr(i,j,klo) = -h_xi*T12 - h_eta*T22 + T23;
732  tau32_arr(i,j,klo) = h_zeta*T23;
733 }
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
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:117
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:104
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:131

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

Referenced by erf_substep_MT(), erf_substep_T(), if(), ERFPhysBCFunct_w::impose_lateral_zvel_bcs(), ERFPhysBCFunct_w::impose_vertical_zvel_bcs(), ParallelFor(), and ERF::project_momenta().

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