ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
SurfaceLayer Class Reference

#include <ERF_SurfaceLayer.H>

Collaboration diagram for SurfaceLayer:

Public Types

enum class  FluxCalcType {
  MOENG = 0 , DONELAN , CUSTOM , BULK_COEFF ,
  ROTATE , RICO
}
 
enum class  ThetaCalcType { ADIABATIC = 0 , HEAT_FLUX , SURFACE_TEMPERATURE }
 
enum class  MoistCalcType { ADIABATIC = 0 , MOISTURE_FLUX , SURFACE_MOISTURE }
 
enum class  RoughCalcType {
  CONSTANT = 0 , CHARNOCK , MODIFIED_CHARNOCK , DONELAN ,
  WAVE_COUPLED
}
 
enum class  PBLHeightCalcType { None , MYNN25 , YSU , MRF }
 

Public Member Functions

 SurfaceLayer (const amrex::Vector< amrex::Geometry > &geom, bool &use_rot_surface_flux, std::string a_pp_prefix, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, const MeshType &a_mesh_type, const TerrainType &a_terrain_type, const TurbChoice &a_turb_choice, amrex::Real start_time, amrex::Real stop_time, amrex::Real bdy_time_interval=0.0)
 
void make_SurfaceLayer_at_level (const int &lev, int nlevs, const amrex::Vector< amrex::MultiFab * > &mfv, std::unique_ptr< amrex::MultiFab > &Theta_prim, std::unique_ptr< amrex::MultiFab > &Qv_prim, std::unique_ptr< amrex::MultiFab > &Qr_prim, std::unique_ptr< amrex::MultiFab > &z_phys_nd, amrex::MultiFab *Hwave, amrex::MultiFab *Lwave, amrex::MultiFab *eddyDiffs, amrex::Vector< amrex::MultiFab * > lsm_data, amrex::Vector< std::string > lsm_data_name, amrex::Vector< amrex::MultiFab * > lsm_flux, amrex::Vector< std::string > lsm_flux_name, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &sst_lev, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &tsk_lev, amrex::Vector< std::unique_ptr< amrex::iMultiFab >> &lmask_lev)
 
void update_fluxes (const int &lev, const amrex::Real &time, amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, const std::unique_ptr< amrex::MultiFab > &walldist, int max_iters=100)
 
template<typename FluxIter >
void compute_fluxes (const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
 
void init_tke_from_ustar (const int &lev, amrex::MultiFab &cons, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, const amrex::Real tkefac=1.0, const amrex::Real zscale=700.0)
 
void impose_SurfaceLayer_bcs (const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Tau_lev, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, const amrex::MultiFab *z_phys)
 
template<typename FluxCalc >
void compute_SurfaceLayer_bcs (const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Tau_lev, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, const amrex::MultiFab *z_phys, const FluxCalc &flux_comp)
 
void fill_tsurf_with_sst_and_tsk (const int &lev, const amrex::Real &time)
 
void fill_qsurf_with_qsat (const int &lev, const amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd)
 
void get_lsm_tsurf (const int &lev)
 
void update_pblh (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const MoistureComponentIndices &moisture_indices)
 
template<typename PBLHeightEstimator >
void compute_pblh (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indice)
 
void read_custom_roughness (const int &lev, const std::string &fname)
 
void update_surf_temp (const amrex::Real &time)
 
void update_mac_ptrs (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
 
amrex::MultiFab * get_u_star (const int &lev)
 
amrex::MultiFab * get_w_star (const int &lev)
 
amrex::MultiFab * get_t_star (const int &lev)
 
amrex::MultiFab * get_q_star (const int &lev)
 
amrex::MultiFab * get_olen (const int &lev)
 
amrex::MultiFab * get_pblh (const int &lev)
 
const amrex::MultiFab * get_mac_avg (const int &lev, int comp)
 
amrex::MultiFab * get_t_surf (const int &lev)
 
void set_t_surf (const int &lev, const amrex::Real tsurf)
 
amrex::MultiFab * get_q_surf (const int &lev)
 
void set_q_surf (const int &lev, const amrex::Real qsurf)
 
amrex::Real get_zref (const int &lev)
 
amrex::MultiFab * get_z0 (const int &lev)
 
bool have_variable_sea_roughness ()
 
amrex::iMultiFab * get_lmask (const int &lev)
 
int lmask_min_reduce (amrex::iMultiFab &lmask, const int &nghost)
 
void update_sst_ptr (const int lev, const int itime, amrex::MultiFab *sst_ptr)
 
void update_tsk_ptr (const int lev, const int itime, amrex::MultiFab *tsk_ptr)
 
template<typename FluxCalc >
void compute_SurfaceLayer_bcs (const int &lev, Vector< const MultiFab * > mfs, Vector< std::unique_ptr< MultiFab >> &Tau_lev, MultiFab *xheat_flux, MultiFab *yheat_flux, MultiFab *zheat_flux, MultiFab *xqv_flux, MultiFab *yqv_flux, MultiFab *zqv_flux, const MultiFab *z_phys, const FluxCalc &flux_comp)
 
template<typename PBLHeightEstimator >
void compute_pblh (const int &lev, Vector< Vector< MultiFab >> &vars, MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indices)
 

Public Attributes

FluxCalcType flux_type {FluxCalcType::MOENG}
 
ThetaCalcType theta_type {ThetaCalcType::ADIABATIC}
 
MoistCalcType moist_type {MoistCalcType::ADIABATIC}
 
RoughCalcType rough_type_land {RoughCalcType::CONSTANT}
 
RoughCalcType rough_type_sea {RoughCalcType::CHARNOCK}
 
PBLHeightCalcType pblh_type {PBLHeightCalcType::None}
 

Private Attributes

amrex::Vector< amrex::Geometry > m_geom
 
bool m_rotate = false
 
amrex::Real m_start_time
 
amrex::Real m_stop_time
 
amrex::Real m_bdy_time_interval
 
bool m_include_wstar = false
 
amrex::Real z0_const {0.1}
 
amrex::Real default_land_surf_temp {300.}
 
amrex::Real surf_temp
 
amrex::Real surf_heating_rate {0}
 
amrex::Real surf_temp_flux {0}
 
amrex::Real default_land_surf_moist {0.}
 
amrex::Real surf_moist
 
amrex::Real surf_moist_flux {0}
 
amrex::Real custom_ustar {0}
 
amrex::Real custom_tstar {0}
 
amrex::Real custom_qstar {0}
 
amrex::Real custom_rhosurf {0}
 
bool specified_rho_surf {false}
 
amrex::Real cnk_a {0.0185}
 
bool cnk_visc {false}
 
amrex::Real depth {30.0}
 
amrex::Vector< amrex::MultiFab > z_0
 
bool m_var_z0 {false}
 
amrex::Real rico_theta_z0 {298.0}
 
amrex::Real rico_qsat_z0 {0.001}
 
bool use_moisture
 
bool m_has_lsm_tsurf = false
 
int m_lsm_tsurf_indx = -1
 
amrex::Real m_Cd = 0.0
 
amrex::Real m_Ch = 0.0
 
amrex::Real m_Cq = 0.0
 
bool m_ignore_sst = false
 
MOSTAverage m_ma
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
 
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
 
amrex::Vector< std::string > m_lsm_data_name
 
amrex::Vector< std::string > m_lsm_flux_name
 
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
 
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
 
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
 
bool m_update_k_rans = false
 
amrex::Real inv_Cmu2 = 0.0
 
amrex::Real theta_ref = 0.0
 

Detailed Description

Abstraction layer for different surface layer schemes (e.g. MOST, Cd)

van der Laan, P., Kelly, M. C., & Sørensen, N. N. (2017). A new k-epsilon model consistent with Monin-Obukhov similarity theory. Wind Energy, 20(3), 479–489. https://doi.org/10.1002/we.2017

Consistent with Dyer (1974) formulation from page 57, Chapter 2, Modeling the vertical ABL structure in Modelling of Atmospheric Flow Fields, Demetri P Lalas and Corrado F Ratto, January 1996, https://doi.org/10.1142/2975.

Member Enumeration Documentation

◆ FluxCalcType

Enumerator
MOENG 

Moeng functional form.

DONELAN 

Donelan functional form.

CUSTOM 

Custom constant flux functional form.

BULK_COEFF 

Bulk transfer coefficient functional form.

ROTATE 

Terrain rotation flux functional form.

RICO 
623  {
624  MOENG = 0, ///< Moeng functional form
625  DONELAN, ///< Donelan functional form
626  CUSTOM, ///< Custom constant flux functional form
627  BULK_COEFF, ///< Bulk transfer coefficient functional form
628  ROTATE, ///< Terrain rotation flux functional form
629  RICO
630  };

◆ MoistCalcType

Enumerator
ADIABATIC 
MOISTURE_FLUX 

Qv-flux specified.

SURFACE_MOISTURE 

Surface Qv specified.

638  {
639  ADIABATIC = 0,
640  MOISTURE_FLUX, ///< Qv-flux specified
641  SURFACE_MOISTURE ///< Surface Qv specified
642  };

◆ PBLHeightCalcType

Enumerator
None 
MYNN25 
YSU 
MRF 
652 { None, MYNN25, YSU, MRF };

◆ RoughCalcType

Enumerator
CONSTANT 

Constant z0.

CHARNOCK 
MODIFIED_CHARNOCK 
DONELAN 
WAVE_COUPLED 
644  {
645  CONSTANT = 0, ///< Constant z0
646  CHARNOCK,
647  MODIFIED_CHARNOCK,
648  DONELAN,
649  WAVE_COUPLED
650  };

◆ ThetaCalcType

Enumerator
ADIABATIC 
HEAT_FLUX 

Heat-flux specified.

SURFACE_TEMPERATURE 

Surface temperature specified.

632  {
633  ADIABATIC = 0,
634  HEAT_FLUX, ///< Heat-flux specified
635  SURFACE_TEMPERATURE ///< Surface temperature specified
636  };

Constructor & Destructor Documentation

◆ SurfaceLayer()

SurfaceLayer::SurfaceLayer ( const amrex::Vector< amrex::Geometry > &  geom,
bool &  use_rot_surface_flux,
std::string  a_pp_prefix,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qv_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  z_phys_nd,
const MeshType &  a_mesh_type,
const TerrainType &  a_terrain_type,
const TurbChoice a_turb_choice,
amrex::Real  start_time,
amrex::Real  stop_time,
amrex::Real  bdy_time_interval = 0.0 
)
inlineexplicit
46  : m_geom(geom),
47  m_rotate(use_rot_surface_flux),
48  m_start_time(start_time),
49  m_stop_time(stop_time),
50  m_bdy_time_interval(bdy_time_interval),
51  m_ma(geom, (z_phys_nd[0] != nullptr), a_pp_prefix, a_mesh_type, a_terrain_type)
52  {
53  // We have a moisture model if Qv_prim is a valid pointer
54  use_moisture = (Qv_prim[0].get());
55 
56  // Get roughness
57  amrex::ParmParse pp("erf");
58  pp.query("most.z0", z0_const);
59 
60  // Specify how to compute the flux
61  if (use_rot_surface_flux) {
63  } else {
64  std::string flux_string_in;
65  std::string flux_string{"moeng"};
66  auto read_flux = pp.query("surface_layer.flux_type", flux_string_in);
67  if (read_flux) {
68  flux_string = amrex::toLower(flux_string_in);
69  }
70  if (flux_string == "donelan") {
72  } else if (flux_string == "moeng") {
74  } else if (flux_string == "rico") {
76  } else if (flux_string == "bulk_coeff") {
78  } else if (flux_string == "custom") {
80  } else {
81  amrex::Abort("Undefined MOST flux type!");
82  }
83  }
84 
85  // Include w* to handle free convection (Beljaars 1995, QJRMS)
86  pp.query("most.include_wstar", m_include_wstar);
87 
88  std::string pblh_string_in;
89  std::string pblh_string{"none"};
90  auto read_pblh = pp.query("most.pblh_calc", pblh_string_in);
91  if (read_pblh) {
92  pblh_string = amrex::toLower(pblh_string_in);
93  }
94  if (pblh_string == "none") {
96  } else if (pblh_string == "mynn25") {
98  } else if (pblh_string == "mynnedmf") {
100  } else if (pblh_string == "ysu") {
102  } else if (pblh_string == "mrf") {
104  } else {
105  amrex::Abort("Undefined PBLH calc type!");
106  }
107 
108  // Get surface temperature
109  auto erf_st = pp.query("most.surf_temp", surf_temp);
110  if (erf_st) { default_land_surf_temp = surf_temp; }
111 
112  // Get surface moisture
113  bool erf_sq = false;
114  if (use_moisture) { erf_sq = pp.query("most.surf_moist", surf_moist); }
115  if (erf_sq) { default_land_surf_moist = surf_moist; }
116 
117  // Custom type user must specify the fluxes
122  pp.get("most.ustar", custom_ustar);
123  pp.get("most.tstar", custom_tstar);
124  pp.get("most.qstar", custom_qstar);
125  pp.query("most.rhosurf", custom_rhosurf);
126  if (custom_qstar != 0) {
127  AMREX_ASSERT_WITH_MESSAGE(use_moisture,
128  "Specified custom MOST qv flux without moisture model!");
129  }
130  amrex::Print() << "Using specified ustar, tstar, qstar for MOST = "
131  << custom_ustar << " " << custom_tstar << " "
132  << custom_qstar << std::endl;
133 
134  // Bulk transfer coefficient (must specify coeffs and surface values)
135  } else if (flux_type == FluxCalcType::BULK_COEFF) {
136  pp.get("most.Cd", m_Cd);
137  pp.get("most.Ch", m_Ch);
138  pp.get("most.Cq", m_Cq);
139  pp.get("most.surf_temp", default_land_surf_temp);
140  pp.get("most.surf_moist", default_land_surf_moist);
141  amrex::Print() << "Using specified Cd, Ch, Cq for MOST = "
142  << m_Cd << " " << m_Ch << " "
143  << m_Cq << std::endl;
144 
145  // Specify surface temperature/moisture or surface flux
146  } else {
147  if (erf_st) {
149  pp.query("most.surf_heating_rate", surf_heating_rate); // [K/h]
150  surf_heating_rate = surf_heating_rate / 3600.0; // [K/s]
151  if (pp.query("most.surf_temp_flux", surf_temp_flux)) {
152  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
153  }
154  } else {
155  pp.query("most.surf_temp_flux", surf_temp_flux);
156 
157  if (pp.query("most.surf_heating_rate", surf_heating_rate)) {
158  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
159  }
160  if (std::abs(surf_temp_flux) >
163  } else {
165  }
166  }
167 
168  if (erf_sq) {
170  } else {
171  pp.query("most.surf_moist_flux", surf_moist_flux);
172  if (std::abs(surf_moist_flux) >
175  } else {
177  }
178  }
179  }
180 
182  {
183  pp.query("most.rico.theta_z0", rico_theta_z0);
184  pp.query("most.rico.qsat_z0", rico_qsat_z0);
185  }
186 
187  // Make sure the inputs file doesn't try to use most.roughness_type
188  std::string bogus_input;
189  if (pp.query("most.roughness_type", bogus_input) > 0) {
190  amrex::Abort("most.roughness_type is deprecated; use "
191  "most.roughness_type_land and/or most.roughness_type_sea");
192  }
193 
194  // Specify how to compute the surface flux over land (if there is any)
195  std::string rough_land_string_in;
196  std::string rough_land_string{"constant"};
197  auto read_rough_land =
198  pp.query("most.roughness_type_land", rough_land_string_in);
199  if (read_rough_land) {
200  rough_land_string = amrex::toLower(rough_land_string_in);
201  }
202  if (rough_land_string == "constant") {
204  } else {
205  amrex::Abort("Undefined MOST roughness type for land!");
206  }
207 
208  // Specify how to compute the surface flux over sea (if there is any)
209  std::string rough_sea_string_in;
210  std::string rough_sea_string{"charnock"};
211  auto read_rough_sea = pp.query("most.roughness_type_sea", rough_sea_string_in);
212  if (read_rough_sea) {
213  rough_sea_string = amrex::toLower(rough_sea_string_in);
214  }
215  if (rough_sea_string == "charnock") {
217  pp.query("most.charnock_constant", cnk_a);
218  pp.query("most.charnock_viscosity", cnk_visc);
219  if (cnk_a > 0) {
220  amrex::Print() << "If there is water, Charnock relation with C_a="
221  << cnk_a << (cnk_visc ? " and viscosity" : "")
222  << " will be used" << std::endl;
223  } else {
224  amrex::Print() << "If there is water, Charnock relation with variable "
225  "Charnock parameter (COARE3.0)"
226  << (cnk_visc ? " and viscosity" : "") << " will be used"
227  << std::endl;
228  }
229  } else if (rough_sea_string == "coare3.0") {
231  amrex::Print() << "If there is water, Charnock relation with variable "
232  "Charnock parameter (COARE3.0)"
233  << (cnk_visc ? " and viscosity" : "") << " will be used"
234  << std::endl;
235  cnk_a = -1;
236  } else if (rough_sea_string == "donelan") {
238  } else if (rough_sea_string == "modified_charnock") {
240  pp.query("most.modified_charnock_depth", depth);
241  } else if (rough_sea_string == "wave_coupled") {
243  } else if (rough_sea_string == "constant") {
245  } else {
246  amrex::Abort("Undefined MOST roughness type for sea!");
247  }
248 
249  // use skin temperature instead of sea-surface temperature
250  // (wrfinput data may have lower resolution SST data)
251  pp.query("most.ignore_sst", m_ignore_sst);
252 
253  // If we're using the RANS k model, then we need to update the dirichlet
254  // BC based on the instantaneous u* and θ*; the turbulence modeling
255  // choices can vary per level but for now, assume that if specified then
256  // all levels are using the same RANS model.
257  m_update_k_rans = (a_turb_choice.rans_type == RANSType::kEqn &&
258  a_turb_choice.dirichlet_k == true);
259  if (m_update_k_rans) {
260  inv_Cmu2 = 1.0 / (a_turb_choice.Cmu0 * a_turb_choice.Cmu0);
261  theta_ref = a_turb_choice.theta_ref;
262  }
263 
264  } // constructor
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
ThetaCalcType theta_type
Definition: ERF_SurfaceLayer.H:655
bool m_include_wstar
Definition: ERF_SurfaceLayer.H:669
bool m_rotate
Definition: ERF_SurfaceLayer.H:664
PBLHeightCalcType pblh_type
Definition: ERF_SurfaceLayer.H:659
bool use_moisture
Definition: ERF_SurfaceLayer.H:692
amrex::Real m_Cq
Definition: ERF_SurfaceLayer.H:698
RoughCalcType rough_type_land
Definition: ERF_SurfaceLayer.H:657
amrex::Real z0_const
Definition: ERF_SurfaceLayer.H:670
amrex::Real cnk_a
Definition: ERF_SurfaceLayer.H:683
amrex::Real m_Ch
Definition: ERF_SurfaceLayer.H:697
amrex::Real surf_temp
Definition: ERF_SurfaceLayer.H:672
amrex::Real rico_qsat_z0
Definition: ERF_SurfaceLayer.H:690
bool m_update_k_rans
Definition: ERF_SurfaceLayer.H:723
amrex::Real surf_moist_flux
Definition: ERF_SurfaceLayer.H:677
RoughCalcType rough_type_sea
Definition: ERF_SurfaceLayer.H:658
amrex::Real surf_moist
Definition: ERF_SurfaceLayer.H:676
bool m_ignore_sst
Definition: ERF_SurfaceLayer.H:700
amrex::Real m_bdy_time_interval
Definition: ERF_SurfaceLayer.H:667
amrex::Real m_stop_time
Definition: ERF_SurfaceLayer.H:666
amrex::Real custom_qstar
Definition: ERF_SurfaceLayer.H:680
amrex::Real custom_rhosurf
Definition: ERF_SurfaceLayer.H:681
amrex::Real m_start_time
Definition: ERF_SurfaceLayer.H:665
@ MOENG
Moeng functional form.
@ BULK_COEFF
Bulk transfer coefficient functional form.
@ CUSTOM
Custom constant flux functional form.
@ ROTATE
Terrain rotation flux functional form.
@ DONELAN
Donelan functional form.
@ SURFACE_MOISTURE
Surface Qv specified.
@ MOISTURE_FLUX
Qv-flux specified.
amrex::Real depth
Definition: ERF_SurfaceLayer.H:685
amrex::Real default_land_surf_moist
Definition: ERF_SurfaceLayer.H:675
amrex::Real rico_theta_z0
Definition: ERF_SurfaceLayer.H:689
amrex::Real surf_temp_flux
Definition: ERF_SurfaceLayer.H:674
amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_SurfaceLayer.H:663
amrex::Real theta_ref
Definition: ERF_SurfaceLayer.H:725
amrex::Real custom_tstar
Definition: ERF_SurfaceLayer.H:679
bool cnk_visc
Definition: ERF_SurfaceLayer.H:684
amrex::Real surf_heating_rate
Definition: ERF_SurfaceLayer.H:673
FluxCalcType flux_type
Definition: ERF_SurfaceLayer.H:654
MoistCalcType moist_type
Definition: ERF_SurfaceLayer.H:656
amrex::Real inv_Cmu2
Definition: ERF_SurfaceLayer.H:724
amrex::Real custom_ustar
Definition: ERF_SurfaceLayer.H:678
amrex::Real m_Cd
Definition: ERF_SurfaceLayer.H:696
amrex::Real default_land_surf_temp
Definition: ERF_SurfaceLayer.H:671
@ SURFACE_TEMPERATURE
Surface temperature specified.
@ HEAT_FLUX
Heat-flux specified.
MOSTAverage m_ma
Definition: ERF_SurfaceLayer.H:702
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
RANSType rans_type
Definition: ERF_TurbStruct.H:413
bool dirichlet_k
Definition: ERF_TurbStruct.H:415
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:391
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:402
Here is the call graph for this function:

Member Function Documentation

◆ compute_fluxes()

template<typename FluxIter >
void SurfaceLayer::compute_fluxes ( const int &  lev,
const int &  max_iters,
const FluxIter &  most_flux,
bool  is_land 
)

Function to compute the fluxes (u^star and t^star) for Monin Obukhov similarity theory

Parameters
[in]levCurrent level
[in]max_itersmaximum iterations to use
[in]most_fluxstructure to iteratively compute ustar and tstar
226 {
227  // Pointers to the computed averages
228  const auto *const tm_ptr = m_ma.get_average(lev,2); // potential temperature
229  const auto *const qvm_ptr = m_ma.get_average(lev,3); // water vapor mixing ratio
230  const auto *const tvm_ptr = m_ma.get_average(lev,4); // virtual potential temperature
231  const auto *const umm_ptr = m_ma.get_average(lev,5); // horizontal velocity magnitude
232  const auto *const zref_ptr = m_ma.get_zref(lev); // reference height
233 
234  for (MFIter mfi(*u_star[lev]); mfi.isValid(); ++mfi)
235  {
236  Box gtbx = mfi.growntilebox();
237 
238  auto u_star_arr = u_star[lev]->array(mfi);
239  auto t_star_arr = t_star[lev]->array(mfi);
240  auto q_star_arr = q_star[lev]->array(mfi);
241  auto t_surf_arr = t_surf[lev]->array(mfi);
242  auto q_surf_arr = q_surf[lev]->array(mfi);
243  auto olen_arr = olen[lev]->array(mfi);
244 
245  const auto tm_arr = tm_ptr->array(mfi);
246  const auto tvm_arr = tvm_ptr->array(mfi);
247  const auto qvm_arr = qvm_ptr->array(mfi);
248  const auto umm_arr = umm_ptr->array(mfi);
249  const auto zref_arr = zref_ptr->array(mfi);
250  const auto z0_arr = z_0[lev].array(mfi);
251 
252  // PBL height if we need to calculate wstar for the Beljaars correction
253  // TODO: can/should we apply this in LES mode?
254  const auto w_star_arr = (m_include_wstar) ? w_star[lev].get()->array(mfi) : Array4<Real> {};
255  const auto pblh_arr = (m_include_wstar) ? pblh[lev].get()->array(mfi) : Array4<Real> {};
256 
257  // Wave properties if they exist
258  const auto Hwave_arr = (m_Hwave_lev[lev]) ? m_Hwave_lev[lev]->array(mfi) : Array4<Real> {};
259  const auto Lwave_arr = (m_Lwave_lev[lev]) ? m_Lwave_lev[lev]->array(mfi) : Array4<Real> {};
260  const auto eta_arr = (m_eddyDiffs_lev[lev]) ? m_eddyDiffs_lev[lev]->array(mfi) : Array4<Real> {};
261 
262  // Land mask array if it exists
263  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
264  Array4<int> {};
265 
266  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
267  {
268  if (( is_land && lmask_arr(i,j,k) == 1) ||
269  (!is_land && lmask_arr(i,j,k) == 0))
270  {
271  most_flux.iterate_flux(i, j, k, max_iters,
272  zref_arr, // set in most average
273  z0_arr, // updated if(!is_land)
274  umm_arr, tm_arr, tvm_arr, qvm_arr,
275  u_star_arr, // updated
276  w_star_arr, // updated if(m_include_wstar)
277  t_star_arr, q_star_arr, // updated
278  t_surf_arr, q_surf_arr, olen_arr, // updated
279  pblh_arr, // updated if(m_include_wstar)
280  Hwave_arr, Lwave_arr, eta_arr);
281  }
282  });
283  }
284 }
amrex::MultiFab * get_zref(const int &lev) const
Definition: ERF_MOSTAverage.H:105
const amrex::MultiFab * get_average(const int &lev, const int &comp) const
Definition: ERF_MOSTAverage.H:102
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
Definition: ERF_SurfaceLayer.H:714
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
Definition: ERF_SurfaceLayer.H:709
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
Definition: ERF_SurfaceLayer.H:706
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
Definition: ERF_SurfaceLayer.H:720
amrex::Vector< amrex::MultiFab > z_0
Definition: ERF_SurfaceLayer.H:686
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
Definition: ERF_SurfaceLayer.H:704
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
Definition: ERF_SurfaceLayer.H:703
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
Definition: ERF_SurfaceLayer.H:710
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
Definition: ERF_SurfaceLayer.H:719
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
Definition: ERF_SurfaceLayer.H:705
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
Definition: ERF_SurfaceLayer.H:721
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
Definition: ERF_SurfaceLayer.H:707
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
Definition: ERF_SurfaceLayer.H:708

◆ compute_pblh() [1/2]

template<typename PBLHeightEstimator >
void SurfaceLayer::compute_pblh ( const int &  lev,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars,
amrex::MultiFab *  z_phys_cc,
const PBLHeightEstimator &  est,
const MoistureComponentIndices moisture_indice 
)

◆ compute_pblh() [2/2]

template<typename PBLHeightEstimator >
void SurfaceLayer::compute_pblh ( const int &  lev,
Vector< Vector< MultiFab >> &  vars,
MultiFab *  z_phys_cc,
const PBLHeightEstimator &  est,
const MoistureComponentIndices moisture_indices 
)
755 {
756  est.compute_pblh(m_geom[lev],z_phys_cc, pblh[lev].get(),
757  vars[lev][Vars::cons],m_lmask_lev[lev][0],
758  moisture_indices);
759 }
@ cons
Definition: ERF_IndexDefines.H:140

◆ compute_SurfaceLayer_bcs() [1/2]

template<typename FluxCalc >
void SurfaceLayer::compute_SurfaceLayer_bcs ( const int &  lev,
amrex::Vector< const amrex::MultiFab * >  mfs,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Tau_lev,
amrex::MultiFab *  xheat_flux,
amrex::MultiFab *  yheat_flux,
amrex::MultiFab *  zheat_flux,
amrex::MultiFab *  xqv_flux,
amrex::MultiFab *  yqv_flux,
amrex::MultiFab *  zqv_flux,
const amrex::MultiFab *  z_phys,
const FluxCalc &  flux_comp 
)

◆ compute_SurfaceLayer_bcs() [2/2]

template<typename FluxCalc >
void SurfaceLayer::compute_SurfaceLayer_bcs ( const int &  lev,
Vector< const MultiFab * >  mfs,
Vector< std::unique_ptr< MultiFab >> &  Tau_lev,
MultiFab *  xheat_flux,
MultiFab *  yheat_flux,
MultiFab *  zheat_flux,
MultiFab *  xqv_flux,
MultiFab *  yqv_flux,
MultiFab *  zqv_flux,
const MultiFab *  z_phys,
const FluxCalc &  flux_comp 
)

Function to calculate MOST fluxes for populating ghost cells.

Parameters
[in]levCurrent level
[in,out]mfsMultifabs to populate
[in]eddyDiffsDiffusion coefficients from turbulence model
[in]flux_compstructure to compute fluxes
369 {
370  bool rotate = m_rotate;
371  const int klo = m_geom[lev].Domain().smallEnd(2);
372  const auto& dxInv = m_geom[lev].InvCellSizeArray();
373  for (MFIter mfi(*mfs[0]); mfi.isValid(); ++mfi)
374  {
375  // Get field arrays
376  const auto cons_arr = mfs[Vars::cons]->array(mfi);
377  const auto velx_arr = mfs[Vars::xvel]->array(mfi);
378  const auto vely_arr = mfs[Vars::yvel]->array(mfi);
379  const auto velz_arr = mfs[Vars::zvel]->array(mfi);
380 
381  // Diffusive stress vars
382  auto t13_arr = Tau_lev[TauType::tau13]->array(mfi);
383  auto t31_arr = (Tau_lev[TauType::tau31]) ? Tau_lev[TauType::tau31]->array(mfi) : Array4<Real>{};
384 
385  auto t23_arr = Tau_lev[TauType::tau23]->array(mfi);
386  auto t32_arr = (Tau_lev[TauType::tau32]) ? Tau_lev[TauType::tau32]->array(mfi) : Array4<Real>{};
387 
388  auto hfx3_arr = zheat_flux->array(mfi);
389  auto qfx3_arr = (zqv_flux) ? zqv_flux->array(mfi) : Array4<Real>{};
390 
391  // Rotated stress vars
392  auto t11_arr = (m_rotate) ? Tau_lev[TauType::tau11]->array(mfi) : Array4<Real>{};
393  auto t22_arr = (m_rotate) ? Tau_lev[TauType::tau22]->array(mfi) : Array4<Real>{};
394  auto t33_arr = (m_rotate) ? Tau_lev[TauType::tau33]->array(mfi) : Array4<Real>{};
395  auto t12_arr = (m_rotate) ? Tau_lev[TauType::tau12]->array(mfi) : Array4<Real>{};
396  auto t21_arr = (m_rotate) ? Tau_lev[TauType::tau21]->array(mfi) : Array4<Real>{};
397 
398  auto hfx1_arr = (m_rotate) ? xheat_flux->array(mfi) : Array4<Real>{};
399  auto hfx2_arr = (m_rotate) ? yheat_flux->array(mfi) : Array4<Real>{};
400  auto qfx1_arr = (m_rotate && xqv_flux) ? xqv_flux->array(mfi) : Array4<Real>{};
401  auto qfx2_arr = (m_rotate && yqv_flux) ? yqv_flux->array(mfi) : Array4<Real>{};
402 
403  // Terrain
404  const auto zphys_arr = (z_phys) ? z_phys->const_array(mfi) : Array4<const Real>{};
405 
406  // Get average arrays
407  const auto *const u_mean = m_ma.get_average(lev,0);
408  const auto *const v_mean = m_ma.get_average(lev,1);
409  const auto *const t_mean = m_ma.get_average(lev,2);
410  const auto *const q_mean = m_ma.get_average(lev,3);
411  const auto *const u_mag_mean = m_ma.get_average(lev,5);
412 
413  const auto um_arr = u_mean->array(mfi);
414  const auto vm_arr = v_mean->array(mfi);
415  const auto tm_arr = t_mean->array(mfi);
416  const auto qm_arr = q_mean->array(mfi);
417  const auto umm_arr = u_mag_mean->array(mfi);
418 
419  // Get derived arrays
420  const auto u_star_arr = u_star[lev]->array(mfi);
421  const auto t_star_arr = t_star[lev]->array(mfi);
422  const auto q_star_arr = q_star[lev]->array(mfi);
423  const auto t_surf_arr = t_surf[lev]->array(mfi);
424  const auto q_surf_arr = q_surf[lev]->array(mfi);
425 
426  // Get LSM fluxes
427  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
428  Array4<int> {};
429  auto lsm_t_flux_arr = Array4<Real> {};
430  auto lsm_q_flux_arr = Array4<Real> {};
431  auto lsm_tau13_arr = Array4<Real> {};
432  auto lsm_tau23_arr = Array4<Real> {};
433  for (int n(0); n<m_lsm_flux_lev[lev].size(); ++n) {
434  if (toLower(m_lsm_flux_name[n]) == "t_flux") { lsm_t_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
435  if (toLower(m_lsm_flux_name[n]) == "q_flux") { lsm_q_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
436  if (toLower(m_lsm_flux_name[n]) == "tau13") { lsm_tau13_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
437  if (toLower(m_lsm_flux_name[n]) == "tau23") { lsm_tau23_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
438  }
439 
440 
441  // Rho*Theta flux
442  //============================================================================
443  Box bx = mfi.tilebox();
444  if (bx.smallEnd(2) != klo) { continue; }
445  bx.makeSlab(2,klo);
446  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
447  {
448  // Valid theta flux from LSM and over land
449  Real Tflux;
450  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
451  if (lsm_t_flux_arr && is_land) {
452  Tflux = lsm_t_flux_arr(i,j,k);
453  } else {
454  Tflux = flux_comp.compute_t_flux(i, j, k,
455  cons_arr, velx_arr, vely_arr,
456  umm_arr, tm_arr, u_star_arr,
457  t_star_arr, t_surf_arr);
458 
459  }
460 
461  // Do scalar flux rotations?
462  if (rotate) {
463  rotate_scalar_flux(i, j, k, Tflux, dxInv, zphys_arr,
464  hfx1_arr, hfx2_arr, hfx3_arr);
465  } else {
466  hfx3_arr(i,j,klo) = Tflux;
467  }
468  });
469 
470  // Rho*Qv flux
471  //============================================================================
472  if (use_moisture) {
473  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
474  {
475  // Valid qv flux from LSM and over land
476  Real Qflux;
477  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
478  if (lsm_q_flux_arr && is_land) {
479  Qflux = lsm_q_flux_arr(i,j,k);
480  } else {
481  Qflux = flux_comp.compute_q_flux(i, j, k,
482  cons_arr, velx_arr, vely_arr,
483  umm_arr, qm_arr, u_star_arr,
484  q_star_arr, q_surf_arr);
485  }
486 
487  // Do scalar flux rotations?
488  if (rotate) {
489  rotate_scalar_flux(i, j, k, Qflux, dxInv, zphys_arr,
490  qfx1_arr, qfx2_arr, qfx3_arr);
491  } else {
492  qfx3_arr(i,j,k) = Qflux;
493  }
494  });
495  } // custom
496 
497  if (!rotate) {
498  // Rho*u flux
499  //============================================================================
500  Box bxx = surroundingNodes(bx,0);
501  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
502  {
503  // Valid tau13 from LSM and over land
504  Real stressx;
505  int is_land_hi = (lmask_arr) ? lmask_arr(i ,j,klo) : 1;
506  int is_land_lo = (lmask_arr) ? lmask_arr(i-1,j,klo) : 1;
507  if (lsm_tau13_arr && (is_land_hi || is_land_lo)) {
508  stressx = 0.;
509  if (!is_land_hi || !is_land_lo) {
510  stressx += 0.5 * flux_comp.compute_u_flux(i, j, k,
511  cons_arr, velx_arr, vely_arr,
512  umm_arr, um_arr, u_star_arr);
513  }
514  if (is_land_hi) {
515  stressx += 0.5 * lsm_tau13_arr(i ,j,k);
516  }
517  if (is_land_lo) {
518  stressx += 0.5 * lsm_tau13_arr(i-1,j,k);
519  }
520  } else {
521  stressx = flux_comp.compute_u_flux(i, j, k,
522  cons_arr, velx_arr, vely_arr,
523  umm_arr, um_arr, u_star_arr);
524  }
525 
526  t13_arr(i,j,k) = stressx;
527  if (t31_arr) { t31_arr(i,j,k) = stressx; }
528  });
529 
530  // Rho*v flux
531  //============================================================================
532  Box bxy = surroundingNodes(bx,1);
533  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
534  {
535  // Valid tau13 from LSM and over land
536  Real stressy;
537  int is_land_hi = (lmask_arr) ? lmask_arr(i,j ,klo) : 1;
538  int is_land_lo = (lmask_arr) ? lmask_arr(i,j-1,klo) : 1;
539  if (lsm_tau23_arr && (is_land_hi || is_land_lo)) {
540  stressy = 0.;
541  if (!is_land_hi || !is_land_lo) {
542  stressy += 0.5 * flux_comp.compute_v_flux(i, j, k,
543  cons_arr, velx_arr, vely_arr,
544  umm_arr, vm_arr, u_star_arr);
545  }
546  if (is_land_hi) {
547  stressy += 0.5 * lsm_tau23_arr(i,j ,k);
548  }
549  if (is_land_lo) {
550  stressy += 0.5 * lsm_tau23_arr(i,j-1,k);
551  }
552  } else {
553  stressy = flux_comp.compute_v_flux(i, j, k,
554  cons_arr, velx_arr, vely_arr,
555  umm_arr, vm_arr, u_star_arr);
556  }
557 
558  t23_arr(i,j,k) = stressy;
559  if (t32_arr) { t32_arr(i,j,k) = stressy; }
560  });
561  } else {
562  // All fluxes with rotation
563  //============================================================================
564  Box bxxy = convert(bx, IntVect(1,1,0));
565  ParallelFor(bxxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
566  {
567  Real stresst = flux_comp.compute_u_flux(i, j, k,
568  cons_arr, velx_arr, vely_arr,
569  umm_arr, um_arr, u_star_arr);
570  rotate_stress_tensor(i, j, k, stresst, dxInv, zphys_arr,
571  velx_arr, vely_arr, velz_arr,
572  t11_arr, t22_arr, t33_arr,
573  t12_arr, t21_arr,
574  t13_arr, t31_arr,
575  t23_arr, t32_arr);
576  });
577  }
578  } // mfiter
579 }
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau32
Definition: ERF_DataStruct.H:32
@ tau31
Definition: ERF_DataStruct.H:32
@ tau21
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
amrex::Real Real
Definition: ERF_ShocInterface.H:19
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)
Definition: ERF_TerrainMetrics.H:609
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)
Definition: ERF_TerrainMetrics.H:630
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
Definition: ERF_SurfaceLayer.H:716
amrex::Vector< std::string > m_lsm_flux_name
Definition: ERF_SurfaceLayer.H:718
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142
Here is the call graph for this function:

◆ fill_qsurf_with_qsat()

void SurfaceLayer::fill_qsurf_with_qsat ( const int &  lev,
const amrex::MultiFab &  cons_in,
const std::unique_ptr< amrex::MultiFab > &  z_phys_nd 
)
666 {
667  // NOTE: We have already tested a moisture model exists
668 
669  // Populate q_surf with qsat over water
670  auto dz = m_geom[lev].CellSize(2);
671  for (MFIter mfi(*q_surf[lev]); mfi.isValid(); ++mfi)
672  {
673  Box gtbx = mfi.growntilebox();
674 
675  auto t_surf_arr = t_surf[lev]->array(mfi);
676  auto q_surf_arr = q_surf[lev]->array(mfi);
677  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
678  Array4<int> {};
679  const auto cons_arr = cons_in.const_array(mfi);
680  const auto z_arr = (z_phys_nd) ? z_phys_nd->const_array(mfi) :
681  Array4<const Real> {};
682 
683  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
684  {
685  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
686  if (!is_land) {
687  auto deltaZ = (z_arr) ? Compute_Zrel_AtCellCenter(i,j,k,z_arr) :
688  0.5*dz;
689  auto Rho = cons_arr(i,j,k,Rho_comp);
690  auto RTh = cons_arr(i,j,k,RhoTheta_comp);
691  auto Qv = cons_arr(i,j,k,RhoQ1_comp) / Rho;
692  auto P_cc = getPgivenRTh(RTh, Qv);
693  P_cc += Rho*CONST_GRAV*deltaZ;
694  P_cc *= 0.01;
695  erf_qsatw(t_surf_arr(i,j,k), P_cc, q_surf_arr(i,j,k));
696  }
697  });
698  }
699  q_surf[lev]->FillBoundary(m_geom[lev].periodicity());
700 }
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
#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_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:166
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)
Definition: ERF_TerrainMetrics.H:387
Here is the call graph for this function:

◆ fill_tsurf_with_sst_and_tsk()

void SurfaceLayer::fill_tsurf_with_sst_and_tsk ( const int &  lev,
const amrex::Real time 
)
584 {
585  int n_times_in_sst = m_sst_lev[lev].size();
586 
587  int n_time_lo, n_time_hi;
588  Real alpha;
589 
590  if (n_times_in_sst > 1) {
591  // Time interpolation
593  int n_time = static_cast<int>( elapsed_time / dT);
594  n_time_lo = n_time;
595  n_time_hi = n_time+1;
596  alpha = (elapsed_time - n_time * dT) / dT;
597  if ((elapsed_time == m_stop_time-m_start_time) && (alpha==0)) {
598  // stop time coincides with final lowinp slice -- don't try to
599  // interpolate from the following time slice
600  n_time -= 1;
601  n_time_lo -= 1;
602  n_time_hi -= 1;
603  alpha = 1.0;
604  }
605  AMREX_ALWAYS_ASSERT( (n_time >= 0) && (n_time < (m_sst_lev[lev].size()-1)));
606  } else {
607  n_time_lo = 0;
608  n_time_hi = 0;
609  alpha = 1.0;
610  }
611  AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
612 
613  Real oma = 1.0 - alpha;
614 
615  // Define a default land surface temperature if we don't read in tsk
617 
618  bool use_tsk = (m_tsk_lev[lev][0]);
619  bool ignore_sst = m_ignore_sst;
620 
621  // Populate t_surf
622  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
623  {
624  Box gtbx = mfi.growntilebox();
625 
626  auto t_surf_arr = t_surf[lev]->array(mfi);
627 
628  const auto sst_lo_arr = m_sst_lev[lev][n_time_lo]->const_array(mfi);
629  const auto sst_hi_arr = m_sst_lev[lev][n_time_hi]->const_array(mfi);
630 
631  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
632  Array4<int> {};
633 
634  if (use_tsk) {
635  const auto tsk_arr = m_tsk_lev[lev][n_time_lo]->const_array(mfi);
636  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
637  {
638  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
639  if (!is_land && !ignore_sst) {
640  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
641  + alpha * sst_hi_arr(i,j,k);
642  } else {
643  t_surf_arr(i,j,k) = tsk_arr(i,j,k);
644  }
645  });
646  } else {
647  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
648  {
649  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
650  if (!is_land) {
651  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
652  + alpha * sst_hi_arr(i,j,k);
653  } else {
654  t_surf_arr(i,j,k) = lst;
655  }
656  });
657  }
658  }
659  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
660 }
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
Definition: ERF_SurfaceLayer.H:712
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
Definition: ERF_SurfaceLayer.H:713

◆ get_lmask()

amrex::iMultiFab* SurfaceLayer::get_lmask ( const int &  lev)
inline
593 { return m_lmask_lev[lev][0]; }

◆ get_lsm_tsurf()

void SurfaceLayer::get_lsm_tsurf ( const int &  lev)
704 {
705  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
706  {
707  Box gtbx = mfi.growntilebox();
708 
709  // NOTE: LSM does not carry lateral ghost cells.
710  // This copies the valid box into the ghost cells.
711  // Fillboundary is called after this to pick up the
712  // interior ghost and periodic directions.
713  Box vbx = mfi.validbox();
714  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
715  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
716 
717  auto t_surf_arr = t_surf[lev]->array(mfi);
718  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
719  Array4<int> {};
720  const auto lsm_arr = m_lsm_data_lev[lev][m_lsm_tsurf_indx]->const_array(mfi);
721 
722  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
723  {
724  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
725  if (is_land) {
726  int li = amrex::min(amrex::max(i, i_lo), i_hi);
727  int lj = amrex::min(amrex::max(j, j_lo), j_hi);
728  t_surf_arr(i,j,k) = lsm_arr(li,lj,k);
729  }
730  });
731  }
732 }
int m_lsm_tsurf_indx
Definition: ERF_SurfaceLayer.H:694
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
Definition: ERF_SurfaceLayer.H:715

◆ get_mac_avg()

const amrex::MultiFab* SurfaceLayer::get_mac_avg ( const int &  lev,
int  comp 
)
inline
577  {
578  return m_ma.get_average(lev, comp);
579  }
Here is the call graph for this function:

◆ get_olen()

amrex::MultiFab* SurfaceLayer::get_olen ( const int &  lev)
inline
572 { return olen[lev].get(); }

◆ get_pblh()

amrex::MultiFab* SurfaceLayer::get_pblh ( const int &  lev)
inline
574 { return pblh[lev].get(); }

◆ get_q_star()

amrex::MultiFab* SurfaceLayer::get_q_star ( const int &  lev)
inline
570 { return q_star[lev].get(); }

◆ get_q_surf()

amrex::MultiFab* SurfaceLayer::get_q_surf ( const int &  lev)
inline
584 { return q_surf[lev].get(); }

◆ get_t_star()

amrex::MultiFab* SurfaceLayer::get_t_star ( const int &  lev)
inline
568 { return t_star[lev].get(); }

◆ get_t_surf()

amrex::MultiFab* SurfaceLayer::get_t_surf ( const int &  lev)
inline
581 { return t_surf[lev].get(); }

◆ get_u_star()

amrex::MultiFab* SurfaceLayer::get_u_star ( const int &  lev)
inline
564 { return u_star[lev].get(); }

◆ get_w_star()

amrex::MultiFab* SurfaceLayer::get_w_star ( const int &  lev)
inline
566 { return w_star[lev].get(); }

◆ get_z0()

amrex::MultiFab* SurfaceLayer::get_z0 ( const int &  lev)
inline
589 { return &z_0[lev]; }

◆ get_zref()

amrex::Real SurfaceLayer::get_zref ( const int &  lev)
inline
587 { return (m_ma.get_zref(lev))->min(0); }
Here is the call graph for this function:

◆ have_variable_sea_roughness()

bool SurfaceLayer::have_variable_sea_roughness ( )
inline
591 { return m_var_z0; }
bool m_var_z0
Definition: ERF_SurfaceLayer.H:687

◆ impose_SurfaceLayer_bcs()

void SurfaceLayer::impose_SurfaceLayer_bcs ( const int &  lev,
amrex::Vector< const amrex::MultiFab * >  mfs,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Tau_lev,
amrex::MultiFab *  xheat_flux,
amrex::MultiFab *  yheat_flux,
amrex::MultiFab *  zheat_flux,
amrex::MultiFab *  xqv_flux,
amrex::MultiFab *  yqv_flux,
amrex::MultiFab *  zqv_flux,
const amrex::MultiFab *  z_phys 
)

Wrapper to impose Monin Obukhov similarity theory fluxes by populating ghost cells.

Parameters
[in]levCurrent level
[in,out]mfsMultifabs to populate
[in]eddyDiffsDiffusion coefficients from turbulence model
305 {
307  moeng_flux flux_comp;
308  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
309  xheat_flux, yheat_flux, zheat_flux,
310  xqv_flux, yqv_flux, zqv_flux,
311  z_phys, flux_comp);
312  } else if (flux_type == FluxCalcType::DONELAN) {
313  donelan_flux flux_comp;
314  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
315  xheat_flux, yheat_flux, zheat_flux,
316  xqv_flux, yqv_flux, zqv_flux,
317  z_phys, flux_comp);
318  } else if (flux_type == FluxCalcType::ROTATE) {
319  rotate_flux flux_comp;
320  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
321  xheat_flux, yheat_flux, zheat_flux,
322  xqv_flux, yqv_flux, zqv_flux,
323  z_phys, flux_comp);
324  } else if (flux_type == FluxCalcType::RICO) {
326  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
327  xheat_flux, yheat_flux, zheat_flux,
328  xqv_flux, yqv_flux, zqv_flux,
329  z_phys, flux_comp);
330  } else if (flux_type == FluxCalcType::BULK_COEFF) {
331  bulk_coeff_flux flux_comp(m_Cd, m_Ch, m_Cq);
332  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
333  xheat_flux, yheat_flux, zheat_flux,
334  xqv_flux, yqv_flux, zqv_flux,
335  z_phys, flux_comp);
336  } else if (flux_type == FluxCalcType::CUSTOM) {
337  custom_flux flux_comp(specified_rho_surf);
338  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
339  xheat_flux, yheat_flux, zheat_flux,
340  xqv_flux, yqv_flux, zqv_flux,
341  z_phys, flux_comp);
342  } else {
343  amrex::Abort("Unknown surface layer flux calculation type");
344  }
345 }
bool specified_rho_surf
Definition: ERF_SurfaceLayer.H:682
void compute_SurfaceLayer_bcs(const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Tau_lev, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, const amrex::MultiFab *z_phys, const FluxCalc &flux_comp)
Definition: ERF_MOSTStress.H:2180
Definition: ERF_MOSTStress.H:2063
Definition: ERF_MOSTStress.H:1931
Definition: ERF_MOSTStress.H:1768
Definition: ERF_MOSTStress.H:2307
Definition: ERF_MOSTStress.H:2435

◆ init_tke_from_ustar()

void SurfaceLayer::init_tke_from_ustar ( const int &  lev,
amrex::MultiFab &  cons,
const std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const amrex::Real  tkefac = 1.0,
const amrex::Real  zscale = 700.0 
)
767 {
768  Print() << "Initializing TKE from surface layer ustar on level " << lev << std::endl;
769 
770  constexpr Real small = 0.01;
771 
772  for (MFIter mfi(cons); mfi.isValid(); ++mfi)
773  {
774  Box gtbx = mfi.tilebox();
775 
776  auto const& u_star_arr = u_star[lev]->const_array(mfi);
777  auto const& z_phys_arr = z_phys_nd->const_array(mfi);
778 
779  auto const& cons_arr = cons.array(mfi);
780 
781  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
782  {
783  Real rho = cons_arr(i, j, k, Rho_comp);
784  Real ust = u_star_arr(i, j, 0);
785  Real tke0 = tkefac * ust * ust; // surface value
786  Real zagl = Compute_Zrel_AtCellCenter(i, j, k, z_phys_arr);
787 
788  // linearly tapering profile -- following WRF, approximate top of
789  // PBL as ustar * zscale
790  cons_arr(i, j, k, RhoKE_comp) = rho * tke0 * std::max(
791  (ust * zscale - zagl) / (std::max(ust, small) * zscale),
792  small);
793  });
794  }
795 }
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
@ rho
Definition: ERF_Kessler.H:22
Here is the call graph for this function:

◆ lmask_min_reduce()

int SurfaceLayer::lmask_min_reduce ( amrex::iMultiFab &  lmask,
const int &  nghost 
)
inline
597  {
598  int lmask_min = amrex::ReduceMin(lmask, nghost, [=] AMREX_GPU_HOST_DEVICE(
599  amrex::Box const& bx, amrex::Array4<int const> const& lm_arr) -> int
600  {
601  int locmin = std::numeric_limits<int>::max();
602  const auto lo = lbound(bx);
603  const auto hi = ubound(bx);
604  for (int j = lo.y; j <= hi.y; ++j) {
605  for (int i = lo.x; i <= hi.x; ++i) {
606  locmin = std::min(locmin, lm_arr(i, j, 0));
607  }
608  }
609  return locmin;
610  });
611 
612  return lmask_min;
613  }

Referenced by make_SurfaceLayer_at_level().

Here is the caller graph for this function:

◆ make_SurfaceLayer_at_level()

void SurfaceLayer::make_SurfaceLayer_at_level ( const int &  lev,
int  nlevs,
const amrex::Vector< amrex::MultiFab * > &  mfv,
std::unique_ptr< amrex::MultiFab > &  Theta_prim,
std::unique_ptr< amrex::MultiFab > &  Qv_prim,
std::unique_ptr< amrex::MultiFab > &  Qr_prim,
std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
amrex::MultiFab *  Hwave,
amrex::MultiFab *  Lwave,
amrex::MultiFab *  eddyDiffs,
amrex::Vector< amrex::MultiFab * >  lsm_data,
amrex::Vector< std::string >  lsm_data_name,
amrex::Vector< amrex::MultiFab * >  lsm_flux,
amrex::Vector< std::string >  lsm_flux_name,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  sst_lev,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  tsk_lev,
amrex::Vector< std::unique_ptr< amrex::iMultiFab >> &  lmask_lev 
)
inline
283  {
284  // Update MOST Average
286  Theta_prim, Qv_prim, Qr_prim,
287  z_phys_nd);
288 
289  // Get CC vars
290  amrex::MultiFab& mf = *(mfv[0]);
291 
292  amrex::ParmParse pp("erf");
293 
294  // Do we have a time-varying surface roughness that needs to be saved?
295  if (lev == 0) {
296  const int nghost = 0; // ghost cells not included
297  int lmask_min = lmask_min_reduce(*lmask_lev[0].get(), nghost);
298  amrex::ParallelDescriptor::ReduceIntMin(lmask_min);
299 
300  m_var_z0 = (lmask_min < 1) & (rough_type_sea != RoughCalcType::CONSTANT);
301  if (m_var_z0) {
302  std::string rough_sea_string{"charnock"};
303  pp.query("most.roughness_type_sea", rough_sea_string);
304  amrex::Print() << "Variable sea roughness (type " << rough_sea_string
305  << ")" << std::endl;
306  }
307  }
308 
309  if (m_eddyDiffs_lev.size() < lev+1) {
310  m_Hwave_lev.resize(nlevs);
311  m_Lwave_lev.resize(nlevs);
312  m_eddyDiffs_lev.resize(nlevs);
313 
314  m_lsm_data_lev.resize(nlevs);
315  m_lsm_flux_lev.resize(nlevs);
316 
317  m_sst_lev.resize(nlevs);
318  m_tsk_lev.resize(nlevs);
319  m_lmask_lev.resize(nlevs);
320 
321  // Size the MOST params for all levels
322  z_0.resize(nlevs);
323  u_star.resize(nlevs);
324  w_star.resize(nlevs);
325  t_star.resize(nlevs);
326  q_star.resize(nlevs);
327  t_surf.resize(nlevs);
328  q_surf.resize(nlevs);
329  olen.resize(nlevs);
330  pblh.resize(nlevs);
331  }
332 
333  // Get pointers to SST,TSK and LANDMASK data
334  int nt_tot_sst = sst_lev.size();
335  m_sst_lev[lev].resize(nt_tot_sst);
336  for (int nt(0); nt < nt_tot_sst; ++nt) {
337  m_sst_lev[lev][nt] = sst_lev[nt].get();
338  }
339  int nt_tot_tsk = static_cast<int>(tsk_lev.size());
340  m_tsk_lev[lev].resize(nt_tot_tsk);
341  for (int nt(0); nt < nt_tot_tsk; ++nt) {
342  m_tsk_lev[lev][nt] = tsk_lev[nt].get();
343  }
344  int nt_tot_lmask = static_cast<int>(lmask_lev.size());
345  m_lmask_lev[lev].resize(nt_tot_lmask);
346  for (int nt(0); nt < nt_tot_lmask; ++nt) {
347  m_lmask_lev[lev][nt] = lmask_lev[nt].get();
348  }
349 
350  // Get pointers to wave data
351  m_Hwave_lev[lev] = Hwave;
352  m_Lwave_lev[lev] = Lwave;
353  m_eddyDiffs_lev[lev] = eddyDiffs;
354 
355  // Get pointers to LSM data and Fluxes
356  int ndata = static_cast<int>(lsm_data.size());
357  int nflux = static_cast<int>(lsm_flux.size());
358  m_lsm_data_name.resize(ndata);
359  m_lsm_data_lev[lev].resize(ndata);
360  m_lsm_flux_name.resize(nflux);
361  m_lsm_flux_lev[lev].resize(nflux);
362  for (int n(0); n < ndata; ++n) {
363  m_lsm_data_name[n] = lsm_data_name[n];
364  m_lsm_data_lev[lev][n] = lsm_data[n];
365  if (amrex::toLower(lsm_data_name[n]) == "theta") {
366  m_has_lsm_tsurf = true;
367  m_lsm_tsurf_indx = n;
368  }
369  }
370  for (int n(0); n < nflux; ++n) {
371  m_lsm_flux_name[n] = lsm_flux_name[n];
372  m_lsm_flux_lev[lev][n] = lsm_flux[n];
373  }
374 
375  // Check if there is a user-specified roughness file to be read
376  std::string fname;
377  bool read_z0 = false;
378  if ( (flux_type == FluxCalcType::MOENG) ||
380  int count = pp.countval("most.roughness_file_name");
381  if (count > 1) {
382  AMREX_ALWAYS_ASSERT(count >= lev+1);
383  pp.query("most.roughness_file_name", fname, lev);
384  read_z0 = true;
385  } else if (count == 1) {
386  if (lev == 0) {
387  pp.query("most.roughness_file_name", fname);
388  } else {
389  // we will interpolate from the coarsest level
390  fname = "";
391  }
392  read_z0 = true;
393  }
394  // else use z0_const
395  }
396 
397  // Attributes for MFs and FABs
398  //--------------------------------------------------------
399  // Create a 2D ba, dm, & ghost cells
400  amrex::BoxArray ba = mf.boxArray();
401  amrex::BoxList bl2d = ba.boxList();
402  for (auto& b : bl2d) {
403  b.setRange(2, 0);
404  }
405  amrex::BoxArray ba2d(std::move(bl2d));
406  const amrex::DistributionMapping& dm = mf.DistributionMap();
407  const int ncomp = 1;
408  amrex::IntVect ng{1,1,0};
409 
410  // Z0 heights FAB
411  //--------------------------------------------------------
412  z_0[lev].define(ba2d, dm, ncomp, ng);
413  z_0[lev].setVal(z0_const);
414  if (read_z0) {
415  read_custom_roughness(lev, fname);
416  }
417 
418  // 2D MFs for U*, T*, T_surf
419  //--------------------------------------------------------
420  u_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
421  u_star[lev]->setVal(1.E34);
422 
423  w_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
424  w_star[lev]->setVal(1.E34);
425 
426  t_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
427  t_star[lev]->setVal(0.0); // default to neutral
428 
429  q_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
430  q_star[lev]->setVal(0.0); // default to dry
431 
432  olen[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
433  olen[lev]->setVal(1.E34);
434 
435  pblh[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
436  pblh[lev]->setVal(1.E34);
437 
438  t_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
439  t_surf[lev]->setVal(default_land_surf_temp);
440 
441  q_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
442  q_surf[lev]->setVal(default_land_surf_moist);
443 
444  // TODO: Do we want an enum struct for indexing?
445 
446  if (m_sst_lev[lev][0] || m_tsk_lev[lev][0] || m_has_lsm_tsurf) {
447  // Valid SST, TSK or LSM data; t_surf set before computing fluxes (avoids
448  // extended lambda capture) Note that land temp will be set from m_tsk_lev
449  // while sea temp will be set from m_sst_lev
451 
452  // Pathways in fill_tsurf_with_sst_and_tsk
453  bool use_tsk = (m_tsk_lev[lev][0]);
454  bool use_sst = (m_sst_lev[lev][0]);
455  amrex::Print() << "Using MOST with specified surface temperature ";
456  // NOTE: SST from the LOW file populates TSK in update_sst_tsk.
457  // So if we have TSK, it contains everything and has been
458  // sanity checked for valid SST values.
459  if (use_tsk) { m_ignore_sst = true; }
460  if (use_tsk) {
461  amrex::Print() << "(land: TSK, ";
462  } else {
463  amrex::Print() << "(land: T0, ";
464  }
465  if (use_tsk && !use_sst) {
466  amrex::Print() << "sea: TSK)" << std::endl;
467  } else {
468  amrex::Print() << "sea: SST)" << std::endl;
469  AMREX_ALWAYS_ASSERT(m_sst_lev[lev][0]);
470  }
471  }
472  }
void make_MOSTAverage_at_level(const int &lev, const amrex::Vector< amrex::MultiFab * > &vars_old, std::unique_ptr< amrex::MultiFab > &Theta_prim, std::unique_ptr< amrex::MultiFab > &Qv_prim, std::unique_ptr< amrex::MultiFab > &Qr_prim, std::unique_ptr< amrex::MultiFab > &z_phys_nd)
Definition: ERF_MOSTAverage.cpp:85
int lmask_min_reduce(amrex::iMultiFab &lmask, const int &nghost)
Definition: ERF_SurfaceLayer.H:595
amrex::Vector< std::string > m_lsm_data_name
Definition: ERF_SurfaceLayer.H:717
bool m_has_lsm_tsurf
Definition: ERF_SurfaceLayer.H:693
void read_custom_roughness(const int &lev, const std::string &fname)
Definition: ERF_SurfaceLayer.cpp:799
@ ng
Definition: ERF_Morrison.H:48
Here is the call graph for this function:

◆ read_custom_roughness()

void SurfaceLayer::read_custom_roughness ( const int &  lev,
const std::string &  fname 
)
801 {
802  // Read the file if we have it
803  if (!fname.empty()) {
804  // Only the ioproc reads the file
805  Gpu::HostVector<Real> m_x,m_y,m_z0;
806  if (ParallelDescriptor::IOProcessor()) {
807  Print()<<"Reading MOST roughness file at level " << lev << " : " << fname << std::endl;
808  std::ifstream file(fname);
809  Real value1,value2,value3;
810  while(file>>value1>>value2>>value3){
811  m_x.push_back(value1);
812  m_y.push_back(value2);
813  m_z0.push_back(value3);
814  }
815  file.close();
816 
817  AMREX_ALWAYS_ASSERT(m_x.size() == m_y.size());
818  AMREX_ALWAYS_ASSERT(m_x.size() == m_z0.size());
819  }
820 
821  // Broadcast the whole domain to every rank
822  int ioproc = ParallelDescriptor::IOProcessorNumber();
823  int nnode = m_x.size();
824  ParallelDescriptor::Bcast(&nnode, 1, ioproc);
825 
826  if (!ParallelDescriptor::IOProcessor()) {
827  m_x.resize(nnode);
828  m_y.resize(nnode);
829  m_z0.resize(nnode);
830  }
831  ParallelDescriptor::Bcast(m_x.data() , nnode, ioproc);
832  ParallelDescriptor::Bcast(m_y.data() , nnode, ioproc);
833  ParallelDescriptor::Bcast(m_z0.data(), nnode, ioproc);
834 
835  // Copy data to the GPU
836  Gpu::DeviceVector<Real> d_x(nnode),d_y(nnode),d_z0(nnode);
837  Gpu::copy(Gpu::hostToDevice, m_x.begin(), m_x.end(), d_x.begin());
838  Gpu::copy(Gpu::hostToDevice, m_y.begin(), m_y.end(), d_y.begin());
839  Gpu::copy(Gpu::hostToDevice, m_z0.begin(), m_z0.end(), d_z0.begin());
840  Real* xp = d_x.data();
841  Real* yp = d_y.data();
842  Real* z0p = d_z0.data();
843 
844  // Each rank populates it's z_0[lev] MultiFab
845  for (MFIter mfi(z_0[lev]); mfi.isValid(); ++mfi)
846  {
847  Box gtbx = mfi.growntilebox();
848 
849  // Populate z_phys data
850  Real tol = 1.0e-4;
851  auto dx = m_geom[lev].CellSizeArray();
852  auto ProbLoArr = m_geom[lev].ProbLoArray();
853  int ilo = m_geom[lev].Domain().smallEnd(0);
854  int jlo = m_geom[lev].Domain().smallEnd(1);
855  int klo = 0;
856  int ihi = m_geom[lev].Domain().bigEnd(0);
857  int jhi = m_geom[lev].Domain().bigEnd(1);
858 
859  Array4<Real> const& z0_arr = z_0[lev].array(mfi);
860  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
861  {
862  // Clip indices for ghost-cells
863  int ii = amrex::min(amrex::max(i,ilo),ihi);
864  int jj = amrex::min(amrex::max(j,jlo),jhi);
865 
866  // Location of nodes
867  Real x = ProbLoArr[0] + ii * dx[0];
868  Real y = ProbLoArr[1] + jj * dx[1];
869  int inode = ii + jj * (ihi-ilo+2); // stride is Nx+1
870  if (std::sqrt(std::pow(x-xp[inode],2)+std::pow(y-yp[inode],2)) < tol) {
871  z0_arr(i,j,klo) = z0p[inode];
872  } else {
873  // Unexpected list order, do brute force search
874  Real z0loc = 0.0;
875  bool found = false;
876  for (int n=0; n<nnode; ++n) {
877  Real delta=std::sqrt(std::pow(x-xp[n],2)+std::pow(y-yp[n],2));
878  if (delta < tol) {
879  found = true;
880  z0loc = z0p[n];
881  break;
882  }
883  }
884  AMREX_ASSERT_WITH_MESSAGE(found, "Location read from terrain file does not match the grid!");
885  amrex::ignore_unused(found);
886  z0_arr(i,j,klo) = z0loc;
887  }
888  });
889  } // mfi
890  } else {
891  AMREX_ALWAYS_ASSERT(lev > 0);
892 
893  Print()<<"Interpolating MOST roughness at level " << lev << std::endl;
894 
895  // Create a BC mapper that uses FOEXTRAP at domain bndry
896  Vector<int> bc_lo(3,ERFBCType::foextrap);
897  Vector<int> bc_hi(3,ERFBCType::foextrap);
898  Vector<BCRec> bcr; bcr.push_back(BCRec(bc_lo.data(),bc_hi.data()));
899 
900  // Create ref ratio
901  IntVect ratio;
902  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
903  ratio[idim] = m_geom[lev].Domain().length(idim) / m_geom[0].Domain().length(idim);
904  }
905 
906  // Create interp object and interpolate from the coarsest grid
907  MFInterpolater* interp = &mf_cell_cons_interp;
908  interp->interp(z_0[0] , 0,
909  z_0[lev], 0,
910  1, z_0[lev].nGrowVect(),
911  m_geom[0], m_geom[lev],
912  m_geom[lev].Domain(),ratio,
913  bcr, 0);
914  }
915 }
@ m_y
Definition: ERF_DataStruct.H:25
@ m_x
Definition: ERF_DataStruct.H:24
@ foextrap
Definition: ERF_IndexDefines.H:208

Referenced by make_SurfaceLayer_at_level().

Here is the caller graph for this function:

◆ set_q_surf()

void SurfaceLayer::set_q_surf ( const int &  lev,
const amrex::Real  qsurf 
)
inline
585 { q_surf[lev]->setVal(qsurf); }

◆ set_t_surf()

void SurfaceLayer::set_t_surf ( const int &  lev,
const amrex::Real  tsurf 
)
inline
582 { t_surf[lev]->setVal(tsurf); }

◆ update_fluxes()

void SurfaceLayer::update_fluxes ( const int &  lev,
const amrex::Real time,
amrex::MultiFab &  cons_in,
const std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const std::unique_ptr< amrex::MultiFab > &  walldist,
int  max_iters = 100 
)

Wrapper to update ustar and tstar for Monin Obukhov similarity theory.

Parameters
[in]levCurrent level
[in]max_itersmaximum iterations to use
18 {
19  // Update with SST/TSK data if we have a valid pointer
20  if (m_sst_lev[lev][0]) {
21  fill_tsurf_with_sst_and_tsk(lev, time);
22  }
23 
24  // Apply heating rate if needed
26  update_surf_temp(time);
27  }
28 
29  // Update qsurf with qsat over sea
30  if (use_moisture) {
31  fill_qsurf_with_qsat(lev, cons_in, z_phys_nd);
32  }
33 
34  // Update land surface temp if we have a valid pointer
35  if (m_has_lsm_tsurf) { get_lsm_tsurf(lev); }
36 
37  // Fill interior ghost cells
38  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
39 
40  // Compute plane averages for all vars (regardless of flux type)
42 
43 
44  // ***************************************************************
45  // Iterate the fluxes if moeng type
46  // First iterate over land -- the only model for surface roughness
47  // over land is RoughCalcType::CONSTANT
48  // ***************************************************************
51  bool is_land = true;
52  // Do we have a constant flux for moisture over land?
53  bool cons_qflux = ( (moist_type == MoistCalcType::MOISTURE_FLUX) ||
57  surface_flux most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
58  compute_fluxes(lev, max_iters, most_flux, is_land);
59  } else {
60  amrex::Abort("Unknown value for rough_type_land");
61  }
64  surface_temp most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
65  compute_fluxes(lev, max_iters, most_flux, is_land);
66  } else {
67  amrex::Abort("Unknown value for rough_type_land");
68  }
69  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
73  compute_fluxes(lev, max_iters, most_flux, is_land);
74  } else {
75  amrex::Abort("Unknown value for rough_type_land");
76  }
77  } else {
78  amrex::Abort("Unknown value for theta_type");
79  }
80  } // MOENG -- LAND
81 
82  // ***************************************************************
83  // Iterate the fluxes if moeng type
84  // Next iterate over sea -- the models for surface roughness
85  // over sea are CHARNOCK, DONELAN, MODIFIED_CHARNOCK or WAVE_COUPLED
86  // ***************************************************************
89  bool is_land = false;
90  // NOTE: Do not allow default to adiabatic over sea (we have Qvs at surface)
91  // Do we have a constant flux for moisture over sea?
92  bool cons_qflux = (moist_type == MoistCalcType::MOISTURE_FLUX);
96  cnk_a, cnk_visc, cons_qflux);
97  compute_fluxes(lev, max_iters, most_flux, is_land);
100  depth, cons_qflux);
101  compute_fluxes(lev, max_iters, most_flux, is_land);
102  } else if (rough_type_sea == RoughCalcType::DONELAN) {
104  cons_qflux);
105  compute_fluxes(lev, max_iters, most_flux, is_land);
108  cons_qflux);
109  compute_fluxes(lev, max_iters, most_flux, is_land);
110  } else {
111  amrex::Abort("Unknown value for rough_type_sea");
112  }
113 
117  cnk_a, cnk_visc, cons_qflux);
118  compute_fluxes(lev, max_iters, most_flux, is_land);
121  depth, cons_qflux);
122  compute_fluxes(lev, max_iters, most_flux, is_land);
123  } else if (rough_type_sea == RoughCalcType::DONELAN) {
125  cons_qflux);
126  compute_fluxes(lev, max_iters, most_flux, is_land);
129  cons_qflux);
130  compute_fluxes(lev, max_iters, most_flux, is_land);
131  } else {
132  amrex::Abort("Unknown value for rough_type_sea");
133  }
134 
135  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
139  cnk_a, cnk_visc);
140  compute_fluxes(lev, max_iters, most_flux, is_land);
143  depth);
144  compute_fluxes(lev, max_iters, most_flux, is_land);
145  } else if (rough_type_sea == RoughCalcType::DONELAN) {
147  compute_fluxes(lev, max_iters, most_flux, is_land);
150  compute_fluxes(lev, max_iters, most_flux, is_land);
151  } else {
152  amrex::Abort("Unknown value for rough_type_sea");
153  }
154  } else {
155  amrex::Abort("Unknown value for theta_type");
156  }
157 
158  } // MOENG -- SEA
159 
161  if (custom_rhosurf > 0) {
162  specified_rho_surf = true;
163  u_star[lev]->setVal(std::sqrt(custom_rhosurf) * custom_ustar);
164  t_star[lev]->setVal(custom_rhosurf * custom_tstar);
165  q_star[lev]->setVal(custom_rhosurf * custom_qstar);
166  } else {
167  u_star[lev]->setVal(custom_ustar);
168  t_star[lev]->setVal(custom_tstar);
169  q_star[lev]->setVal(custom_qstar);
170  }
171  }
172 
173  if (m_update_k_rans) {
174  const bool use_ref_theta = (theta_ref > 0);
175  const Real l_inv_theta0 = (use_ref_theta) ? 1.0 / theta_ref : 1.0;
176  const Real l_inv_Cmu2 = inv_Cmu2;
177 
178  for (MFIter mfi(*u_star[lev]); mfi.isValid(); ++mfi)
179  {
180  Box gtbx = mfi.growntilebox();
181 
182  auto cons_arr = cons_in.array(mfi);
183  const auto& u_star_arr = u_star[lev]->const_array(mfi);
184  const auto& t_star_arr = t_star[lev]->const_array(mfi);
185  const auto& dist_arr = walldist->const_array(mfi);
186 
187  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
188  {
189  Real rho = cons_arr(i,j,k,Rho_comp);
190  if (t_star_arr(i,j,k) < -1e-8) {
191  // Only destabilizing buoyancy flux affects the boundary k
192  // tstar < 0 ==> B > 0
193  Real B = -CONST_GRAV * l_inv_theta0 * u_star_arr(i,j,k) * t_star_arr(i,j,k);
194  if (!use_ref_theta) {
195  B *= cons_arr(i,j,k,Rho_comp) /
196  cons_arr(i,j,k,RhoTheta_comp);
197  }
198 
199  // Axell & Liungman 2001, Eqn. 16
200  cons_arr(i,j,k,RhoKE_comp) = rho * l_inv_Cmu2 *
201  std::pow(
202  u_star_arr(i,j,k) * u_star_arr(i,j,k) * u_star_arr(i,j,k)
203  + KAPPA * B * dist_arr(i,j,k),
204  2.0/3.0);
205  } else {
206  cons_arr(i,j,k,RhoKE_comp) = rho * l_inv_Cmu2 * u_star_arr(i,j,k) * u_star_arr(i,j,k);
207  }
208  });
209  }
210  }
211 }
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:20
const bool use_ref_theta
Definition: ERF_SetupDiff.H:17
const Real l_inv_theta0
Definition: ERF_SetupDiff.H:18
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:732
void update_surf_temp(const amrex::Real &time)
Definition: ERF_SurfaceLayer.H:543
void get_lsm_tsurf(const int &lev)
Definition: ERF_SurfaceLayer.cpp:703
void fill_qsurf_with_qsat(const int &lev, const amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd)
Definition: ERF_SurfaceLayer.cpp:663
void compute_fluxes(const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
Definition: ERF_SurfaceLayer.cpp:222
void fill_tsurf_with_sst_and_tsk(const int &lev, const amrex::Real &time)
Definition: ERF_SurfaceLayer.cpp:582
Definition: ERF_MOSTStress.H:189
Definition: ERF_MOSTStress.H:359
Definition: ERF_MOSTStress.H:278
Definition: ERF_MOSTStress.H:429
Definition: ERF_MOSTStress.H:140
Definition: ERF_MOSTStress.H:610
Definition: ERF_MOSTStress.H:840
Definition: ERF_MOSTStress.H:729
Definition: ERF_MOSTStress.H:948
Definition: ERF_MOSTStress.H:505
Definition: ERF_MOSTStress.H:1199
Definition: ERF_MOSTStress.H:1523
Definition: ERF_MOSTStress.H:1365
Definition: ERF_MOSTStress.H:1641
Definition: ERF_MOSTStress.H:1065

◆ update_mac_ptrs()

void SurfaceLayer::update_mac_ptrs ( const int &  lev,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars_old,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Theta_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qv_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qr_prim 
)
inline
560  {
561  m_ma.update_field_ptrs(lev, vars_old, Theta_prim, Qv_prim, Qr_prim);
562  }
void update_field_ptrs(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
Definition: ERF_MOSTAverage.cpp:265
Here is the call graph for this function:

◆ update_pblh()

void SurfaceLayer::update_pblh ( const int &  lev,
amrex::Vector< amrex::Vector< amrex::MultiFab >> &  vars,
amrex::MultiFab *  z_phys_cc,
const MoistureComponentIndices moisture_indices 
)
739 {
741  MYNNPBLH estimator;
742  compute_pblh(lev, vars, z_phys_cc, estimator, moisture_indices);
744  amrex::Error("YSU/MRF PBLH calc not implemented yet");
745  }
746 }
void compute_pblh(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indice)
Definition: ERF_PBLHeight.H:8

◆ update_sst_ptr()

void SurfaceLayer::update_sst_ptr ( const int  lev,
const int  itime,
amrex::MultiFab *  sst_ptr 
)
inline
615  {
616  m_sst_lev[lev][itime] = sst_ptr;
617  }

◆ update_surf_temp()

void SurfaceLayer::update_surf_temp ( const amrex::Real time)
inline
544  {
545  if (surf_heating_rate != 0) {
546  int nlevs = static_cast<int>(m_geom.size());
547  for (int lev = 0; lev < nlevs; lev++) {
548  t_surf[lev]->setVal(surf_temp + surf_heating_rate * time);
549  amrex::Print() << "Surface temp at t=" << time << ": "
550  << surf_temp + surf_heating_rate * time << std::endl;
551  }
552  }
553  }

◆ update_tsk_ptr()

void SurfaceLayer::update_tsk_ptr ( const int  lev,
const int  itime,
amrex::MultiFab *  tsk_ptr 
)
inline
619  {
620  m_tsk_lev[lev][itime] = tsk_ptr;
621  }

Member Data Documentation

◆ cnk_a

amrex::Real SurfaceLayer::cnk_a {0.0185}
private

Referenced by SurfaceLayer().

◆ cnk_visc

bool SurfaceLayer::cnk_visc {false}
private

Referenced by SurfaceLayer().

◆ custom_qstar

amrex::Real SurfaceLayer::custom_qstar {0}
private

Referenced by SurfaceLayer().

◆ custom_rhosurf

amrex::Real SurfaceLayer::custom_rhosurf {0}
private

Referenced by SurfaceLayer().

◆ custom_tstar

amrex::Real SurfaceLayer::custom_tstar {0}
private

Referenced by SurfaceLayer().

◆ custom_ustar

amrex::Real SurfaceLayer::custom_ustar {0}
private

Referenced by SurfaceLayer().

◆ default_land_surf_moist

amrex::Real SurfaceLayer::default_land_surf_moist {0.}
private

◆ default_land_surf_temp

amrex::Real SurfaceLayer::default_land_surf_temp {300.}
private

◆ depth

amrex::Real SurfaceLayer::depth {30.0}
private

Referenced by SurfaceLayer().

◆ flux_type

FluxCalcType SurfaceLayer::flux_type {FluxCalcType::MOENG}

◆ inv_Cmu2

amrex::Real SurfaceLayer::inv_Cmu2 = 0.0
private

Referenced by SurfaceLayer().

◆ m_bdy_time_interval

amrex::Real SurfaceLayer::m_bdy_time_interval
private

◆ m_Cd

amrex::Real SurfaceLayer::m_Cd = 0.0
private

Referenced by SurfaceLayer().

◆ m_Ch

amrex::Real SurfaceLayer::m_Ch = 0.0
private

Referenced by SurfaceLayer().

◆ m_Cq

amrex::Real SurfaceLayer::m_Cq = 0.0
private

Referenced by SurfaceLayer().

◆ m_eddyDiffs_lev

amrex::Vector<amrex::MultiFab*> SurfaceLayer::m_eddyDiffs_lev
private

◆ m_geom

amrex::Vector<amrex::Geometry> SurfaceLayer::m_geom
private

Referenced by update_surf_temp().

◆ m_has_lsm_tsurf

bool SurfaceLayer::m_has_lsm_tsurf = false
private

◆ m_Hwave_lev

amrex::Vector<amrex::MultiFab*> SurfaceLayer::m_Hwave_lev
private

◆ m_ignore_sst

bool SurfaceLayer::m_ignore_sst = false
private

◆ m_include_wstar

bool SurfaceLayer::m_include_wstar = false
private

Referenced by SurfaceLayer().

◆ m_lmask_lev

amrex::Vector<amrex::Vector<amrex::iMultiFab*> > SurfaceLayer::m_lmask_lev
private

◆ m_lsm_data_lev

amrex::Vector<amrex::Vector<amrex::MultiFab*> > SurfaceLayer::m_lsm_data_lev
private

◆ m_lsm_data_name

amrex::Vector<std::string> SurfaceLayer::m_lsm_data_name
private

◆ m_lsm_flux_lev

amrex::Vector<amrex::Vector<amrex::MultiFab*> > SurfaceLayer::m_lsm_flux_lev
private

◆ m_lsm_flux_name

amrex::Vector<std::string> SurfaceLayer::m_lsm_flux_name
private

◆ m_lsm_tsurf_indx

int SurfaceLayer::m_lsm_tsurf_indx = -1
private

◆ m_Lwave_lev

amrex::Vector<amrex::MultiFab*> SurfaceLayer::m_Lwave_lev
private

◆ m_ma

MOSTAverage SurfaceLayer::m_ma
private

◆ m_rotate

bool SurfaceLayer::m_rotate = false
private

◆ m_sst_lev

amrex::Vector<amrex::Vector<amrex::MultiFab*> > SurfaceLayer::m_sst_lev
private

◆ m_start_time

amrex::Real SurfaceLayer::m_start_time
private

◆ m_stop_time

amrex::Real SurfaceLayer::m_stop_time
private

◆ m_tsk_lev

amrex::Vector<amrex::Vector<amrex::MultiFab*> > SurfaceLayer::m_tsk_lev
private

◆ m_update_k_rans

bool SurfaceLayer::m_update_k_rans = false
private

Referenced by SurfaceLayer().

◆ m_var_z0

bool SurfaceLayer::m_var_z0 {false}
private

◆ moist_type

MoistCalcType SurfaceLayer::moist_type {MoistCalcType::ADIABATIC}

Referenced by SurfaceLayer().

◆ olen

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::olen
private

◆ pblh

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::pblh
private

◆ pblh_type

PBLHeightCalcType SurfaceLayer::pblh_type {PBLHeightCalcType::None}

Referenced by SurfaceLayer().

◆ q_star

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::q_star
private

◆ q_surf

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::q_surf
private

◆ rico_qsat_z0

amrex::Real SurfaceLayer::rico_qsat_z0 {0.001}
private

Referenced by SurfaceLayer().

◆ rico_theta_z0

amrex::Real SurfaceLayer::rico_theta_z0 {298.0}
private

Referenced by SurfaceLayer().

◆ rough_type_land

RoughCalcType SurfaceLayer::rough_type_land {RoughCalcType::CONSTANT}

Referenced by SurfaceLayer().

◆ rough_type_sea

RoughCalcType SurfaceLayer::rough_type_sea {RoughCalcType::CHARNOCK}

◆ specified_rho_surf

bool SurfaceLayer::specified_rho_surf {false}
private

◆ surf_heating_rate

amrex::Real SurfaceLayer::surf_heating_rate {0}
private

Referenced by SurfaceLayer(), and update_surf_temp().

◆ surf_moist

amrex::Real SurfaceLayer::surf_moist
private

Referenced by SurfaceLayer().

◆ surf_moist_flux

amrex::Real SurfaceLayer::surf_moist_flux {0}
private

Referenced by SurfaceLayer().

◆ surf_temp

amrex::Real SurfaceLayer::surf_temp
private

Referenced by SurfaceLayer(), and update_surf_temp().

◆ surf_temp_flux

amrex::Real SurfaceLayer::surf_temp_flux {0}
private

Referenced by SurfaceLayer().

◆ t_star

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::t_star
private

◆ t_surf

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::t_surf
private

◆ theta_ref

amrex::Real SurfaceLayer::theta_ref = 0.0
private

Referenced by SurfaceLayer().

◆ theta_type

◆ u_star

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::u_star
private

◆ use_moisture

bool SurfaceLayer::use_moisture
private

Referenced by SurfaceLayer().

◆ w_star

amrex::Vector<std::unique_ptr<amrex::MultiFab> > SurfaceLayer::w_star
private

◆ z0_const

amrex::Real SurfaceLayer::z0_const {0.1}
private

◆ z_0

amrex::Vector<amrex::MultiFab> SurfaceLayer::z_0
private

The documentation for this class was generated from the following files: