ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ApplySurfaceTreatment_BulkCoeff.cpp File Reference
#include <AMReX_MultiFab.H>
#include <ERF_SrcHeaders.H>
#include <AMReX_ParmParse.H>
Include dependency graph for ERF_ApplySurfaceTreatment_BulkCoeff.cpp:

Functions

void ApplySurfaceTreatment_BulkCoeff_Mom (const Box &tbx, const Box &tby, const Array4< Real > &rho_u_rhs, const Array4< Real > &rho_v_rhs, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_nd, const Array4< const Real > &surface_state_arr)
 
void ApplySurfaceTreatment_BulkCoeff_CC (const Box &bx, const Array4< Real > &cell_rhs, const Array4< const Real > &cons_state, const Array4< const Real > &z_phys_cc, const Array4< const Real > &surface_state_arr)
 

Function Documentation

◆ ApplySurfaceTreatment_BulkCoeff_CC()

void ApplySurfaceTreatment_BulkCoeff_CC ( const Box &  bx,
const Array4< Real > &  cell_rhs,
const Array4< const Real > &  cons_state,
const Array4< const Real > &  z_phys_cc,
const Array4< const Real > &  surface_state_arr 
)
65 {
66  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
67  if(k == 0) {
68  Real dz = z_phys_cc(i,j,1)-z_phys_cc(i,j,0);
69  Real ls_mask = surface_state_arr(i,j,0);
70  Real Ch = 0.0015*(1.0-ls_mask);
71  Real Ce = 0.0015*(1.0-ls_mask);
72 
73  Real rho = cons_state(i,j,k, Rho_comp);
74  Real rhotheta = cons_state(i,j,k, RhoTheta_comp);
75  Real rhoqv = cons_state(i,j,k, RhoQ1_comp);
76  Real theta = rhotheta/rho;
77  Real qv = rhoqv/rho;
78  Real temp = getTgivenRandRTh(rho, rhotheta, qv);
79  Real uvel = cons_state(i,j,k,1)/cons_state(i,j,k,0);
80  Real vvel = cons_state(i,j,k,2)/cons_state(i,j,k,0);
81  Real velmag = std::sqrt(uvel*uvel + vvel*vvel);
82  Real dT = max(0.0, 301.0 - temp);
83  Real dq = max(0.0, 0.024 - qv);
84  cell_rhs(i, j, k, RhoTheta_comp) += (theta/(1005.0*temp))*rho*Ch*velmag*dT/dz;
85  cell_rhs(i, j, k, RhoQ1_comp) += rho*Ce*velmag*dq/dz;
86  }
87  });
88 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ theta
Definition: ERF_MM5.H:20
@ rho
Definition: ERF_Kessler.H:22
@ qv
Definition: ERF_Kessler.H:28

Referenced by make_sources().

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

◆ ApplySurfaceTreatment_BulkCoeff_Mom()

void ApplySurfaceTreatment_BulkCoeff_Mom ( const Box &  tbx,
const Box &  tby,
const Array4< Real > &  rho_u_rhs,
const Array4< Real > &  rho_v_rhs,
const Array4< const Real > &  rho_u,
const Array4< const Real > &  rho_v,
const Array4< const Real > &  cons_state,
const Array4< const Real > &  z_phys_nd,
const Array4< const Real > &  surface_state_arr 
)
18 {
19  int ndrag = 1;
20  Real Cd_sea = 0.001;
21  Real Cd_land = 0.01;
22 
23  ParallelFor(tbx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
24  {
25  if(k <= ndrag) {
26  Real dz = z_phys_nd(i,j,1) - z_phys_nd(i,j,0);
27  Real ls_mask = (surface_state_arr(i-1,j,0)+surface_state_arr(i,j,0))/2.0;
28  Real weight = 1.0 - Real(k)/ndrag; // linear decrease
29  Real fac = k==0? 1.0 : 0.0;
30  Real Cd = fac*Cd_sea*(1.0-ls_mask) + Cd_land*ls_mask;
31  Real rho_for_u = (cons_state(i-1,j,k,0)+cons_state(i,j,k,0))/2.0;
32  Real rho_for_v = (cons_state(i,j-1,k,0)+cons_state(i,j,k,0))/2.0;
33  Real uvel = rho_u(i,j,k)/rho_for_u;
34  Real vvel = rho_v(i,j,k)/rho_for_v;
35  Real velmag = std::sqrt(uvel*uvel + vvel*vvel);
36  rho_u_rhs(i, j, k) += -1.0*weight*Cd*velmag*rho_for_u*uvel/dz;
37  }
38  });
39 
40 
41  ParallelFor(tby, [=] AMREX_GPU_DEVICE(int i, int j, int k)
42  {
43  if(k <= ndrag) {
44  Real dz = z_phys_nd(i,j,1) - z_phys_nd(i,j,0);
45  Real ls_mask = (surface_state_arr(i,j-1,0)+surface_state_arr(i,j,0))/2.0;
46  Real fac = k==0? 1.0 : 0.0;
47  Real Cd = fac*Cd_sea*(1.0-ls_mask) + Cd_land*ls_mask;
48  Real weight = 1.0 - Real(k)/ndrag; // linear decrease
49  Real rho_for_u = (cons_state(i-1,j,k,0)+cons_state(i,j,k,0))/2.0;
50  Real rho_for_v = (cons_state(i,j-1,k,0)+cons_state(i,j,k,0))/2.0;
51  Real uvel = rho_u(i,j,k)/rho_for_u;
52  Real vvel = rho_v(i,j,k)/rho_for_v;
53  Real velmag = std::sqrt(uvel*uvel + vvel*vvel);
54  rho_v_rhs(i, j, k) += -1.0*weight*Cd*velmag*rho_for_v*vvel/dz;
55  }
56  });
57 }

Referenced by make_mom_sources().

Here is the caller graph for this function: