ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
SurfaceLayer Class Reference

#include <ERF_SurfaceLayer.H>

Collaboration diagram for SurfaceLayer:

Public Types

enum class  FluxCalcType { MOENG = 0 , DONELAN , CUSTOM , ROTATE }
 
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 TerrainType &a_terrain_type, amrex::Real start_bdy_time=0.0, 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< amrex::MultiFab * > lsm_flux, 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, int max_iters=500)
 
template<typename FluxIter >
void compute_fluxes (const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
 
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 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 int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
 
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 int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
 
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)
 
amrex::MultiFab * get_q_surf (const int &lev)
 
amrex::Real get_zref ()
 
amrex::FArrayBox * 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)
 
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, int RhoQv_comp, int RhoQc_comp, int RhoQr_comp)
 

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_bdy_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 cnk_a {0.0185}
 
bool cnk_visc {false}
 
amrex::Real depth {30.0}
 
amrex::Vector< amrex::FArrayBox > z_0
 
bool m_var_z0 {false}
 
bool use_moisture
 
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< amrex::MultiFab * > m_Hwave_lev
 
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
 
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
 

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.

ROTATE 

Terrain rotation flux functional form.

524  {
525  MOENG = 0, ///< Moeng functional form
526  DONELAN, ///< Donelan functional form
527  CUSTOM, ///< Custom constant flux functional form
528  ROTATE ///< Terrain rotation flux functional form
529  };

◆ MoistCalcType

Enumerator
ADIABATIC 
MOISTURE_FLUX 

Qv-flux specified.

SURFACE_MOISTURE 

Surface Qv specified.

537  {
538  ADIABATIC = 0,
539  MOISTURE_FLUX, ///< Qv-flux specified
540  SURFACE_MOISTURE ///< Surface Qv specified
541  };

◆ PBLHeightCalcType

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

◆ RoughCalcType

Enumerator
CONSTANT 

Constant z0.

CHARNOCK 
MODIFIED_CHARNOCK 
DONELAN 
WAVE_COUPLED 
543  {
544  CONSTANT = 0, ///< Constant z0
545  CHARNOCK,
546  MODIFIED_CHARNOCK,
547  DONELAN,
548  WAVE_COUPLED
549  };

◆ ThetaCalcType

Enumerator
ADIABATIC 
HEAT_FLUX 

Heat-flux specified.

SURFACE_TEMPERATURE 

Surface temperature specified.

531  {
532  ADIABATIC = 0,
533  HEAT_FLUX, ///< Heat-flux specified
534  SURFACE_TEMPERATURE ///< Surface temperature specified
535  };

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 TerrainType &  a_terrain_type,
amrex::Real  start_bdy_time = 0.0,
amrex::Real  bdy_time_interval = 0.0 
)
inlineexplicit
42  : m_geom(geom),
43  m_rotate(use_rot_surface_flux),
44  m_start_bdy_time(start_bdy_time),
45  m_bdy_time_interval(bdy_time_interval),
46  m_ma(geom, (z_phys_nd[0] != nullptr), a_pp_prefix, a_terrain_type)
47  {
48  // We have a moisture model if Qv_prim is a valid pointer
49  use_moisture = (Qv_prim[0].get());
50 
51  // Get roughness
52  amrex::ParmParse pp("erf");
53  pp.query("most.z0", z0_const);
54 
55  // Specify how to compute the flux
56  if (use_rot_surface_flux) {
58  } else {
59  std::string flux_string_in;
60  std::string flux_string{"moeng"};
61  auto read_flux = pp.query("surface_layer.flux_type", flux_string_in);
62  if (read_flux) {
63  flux_string = amrex::toLower(flux_string_in);
64  }
65  if (flux_string == "donelan") {
67  } else if (flux_string == "moeng") {
69  } else if (flux_string == "custom") {
71  } else {
72  amrex::Abort("Undefined MOST flux type!");
73  }
74  }
75 
76  // Include w* to handle free convection (Beljaars 1995, QJRMS)
77  pp.query("most.include_wstar", m_include_wstar);
78 
79  std::string pblh_string_in;
80  std::string pblh_string{"none"};
81  auto read_pblh = pp.query("most.pblh_calc", pblh_string_in);
82  if (read_pblh) {
83  pblh_string = amrex::toLower(pblh_string_in);
84  }
85  if (pblh_string == "none") {
87  } else if (pblh_string == "mynn25") {
89  } else if (pblh_string == "mynnedmf") {
91  } else if (pblh_string == "ysu") {
93  } else if (pblh_string == "mrf") {
95  } else {
96  amrex::Abort("Undefined PBLH calc type!");
97  }
98 
99  // Get surface temperature
100  auto erf_st = pp.query("most.surf_temp", surf_temp);
101  if (erf_st) { default_land_surf_temp = surf_temp; }
102 
103  // Get surface moisture
104  bool erf_sq = false;
105  if (use_moisture) { erf_sq = pp.query("most.surf_moist", surf_moist); }
106  if (erf_sq) { default_land_surf_moist = surf_moist; }
107 
108  // Custom type user must specify the fluxes
112  pp.query("most.ustar", custom_ustar);
113  pp.query("most.tstar", custom_tstar);
114  pp.query("most.qstar", custom_qstar);
115  if (custom_qstar != 0) {
116  AMREX_ASSERT_WITH_MESSAGE(use_moisture,
117  "Specified custom MOST qv flux without moisture model!");
118  }
119  amrex::Print() << "Using specified ustar, tstar, qstar for MOST = "
120  << custom_ustar << " " << custom_tstar << " "
121  << custom_qstar << std::endl;
122 
123  // Specify surface temperature/moisture or surface flux
124  } else {
125  if (erf_st) {
127  pp.query("most.surf_heating_rate", surf_heating_rate); // [K/h]
128  surf_heating_rate = surf_heating_rate / 3600.0; // [K/s]
129  if (pp.query("most.surf_temp_flux", surf_temp_flux)) {
130  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
131  }
132  } else {
133  pp.query("most.surf_temp_flux", surf_temp_flux);
134 
135  if (pp.query("most.surf_heating_rate", surf_heating_rate)) {
136  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
137  }
138  if (std::abs(surf_temp_flux) >
141  } else {
143  }
144 
145  }
146 
147  if (erf_sq) {
149  } else {
150  pp.query("most.surf_moist_flux", surf_moist_flux);
151  if (std::abs(surf_moist_flux) >
154  } else {
156  }
157  }
158  }
159 
160  // Make sure the inputs file doesn't try to use most.roughness_type
161  std::string bogus_input;
162  if (pp.query("most.roughness_type", bogus_input) > 0) {
163  amrex::Abort("most.roughness_type is deprecated; use "
164  "most.roughness_type_land and/or most.roughness_type_sea");
165  }
166 
167  // Specify how to compute the surface flux over land (if there is any)
168  std::string rough_land_string_in;
169  std::string rough_land_string{"constant"};
170  auto read_rough_land =
171  pp.query("most.roughness_type_land", rough_land_string_in);
172  if (read_rough_land) {
173  rough_land_string = amrex::toLower(rough_land_string_in);
174  }
175  if (rough_land_string == "constant") {
177  } else {
178  amrex::Abort("Undefined MOST roughness type for land!");
179  }
180 
181  // Specify how to compute the surface flux over sea (if there is any)
182  std::string rough_sea_string_in;
183  std::string rough_sea_string{"charnock"};
184  auto read_rough_sea = pp.query("most.roughness_type_sea", rough_sea_string_in);
185  if (read_rough_sea) {
186  rough_sea_string = amrex::toLower(rough_sea_string_in);
187  }
188  if (rough_sea_string == "charnock") {
190  pp.query("most.charnock_constant", cnk_a);
191  pp.query("most.charnock_viscosity", cnk_visc);
192  if (cnk_a > 0) {
193  amrex::Print() << "If there is water, Charnock relation with C_a="
194  << cnk_a << (cnk_visc ? " and viscosity" : "")
195  << " will be used" << std::endl;
196  } else {
197  amrex::Print() << "If there is water, Charnock relation with variable "
198  "Charnock parameter (COARE3.0)"
199  << (cnk_visc ? " and viscosity" : "") << " will be used"
200  << std::endl;
201  }
202  } else if (rough_sea_string == "coare3.0") {
204  amrex::Print() << "If there is water, Charnock relation with variable "
205  "Charnock parameter (COARE3.0)"
206  << (cnk_visc ? " and viscosity" : "") << " will be used"
207  << std::endl;
208  cnk_a = -1;
209  } else if (rough_sea_string == "donelan") {
211  } else if (rough_sea_string == "modified_charnock") {
213  pp.query("most.modified_charnock_depth", depth);
214  } else if (rough_sea_string == "wave_coupled") {
216  } else if (rough_sea_string == "constant") {
218  } else {
219  amrex::Abort("Undefined MOST roughness type for sea!");
220  }
221  } // constructor
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
ThetaCalcType theta_type
Definition: ERF_SurfaceLayer.H:554
bool m_include_wstar
Definition: ERF_SurfaceLayer.H:567
bool m_rotate
Definition: ERF_SurfaceLayer.H:563
PBLHeightCalcType pblh_type
Definition: ERF_SurfaceLayer.H:558
bool use_moisture
Definition: ERF_SurfaceLayer.H:585
RoughCalcType rough_type_land
Definition: ERF_SurfaceLayer.H:556
amrex::Real z0_const
Definition: ERF_SurfaceLayer.H:568
amrex::Real cnk_a
Definition: ERF_SurfaceLayer.H:579
amrex::Real surf_temp
Definition: ERF_SurfaceLayer.H:570
amrex::Real m_start_bdy_time
Definition: ERF_SurfaceLayer.H:564
amrex::Real surf_moist_flux
Definition: ERF_SurfaceLayer.H:575
RoughCalcType rough_type_sea
Definition: ERF_SurfaceLayer.H:557
amrex::Real surf_moist
Definition: ERF_SurfaceLayer.H:574
amrex::Real m_bdy_time_interval
Definition: ERF_SurfaceLayer.H:565
amrex::Real custom_qstar
Definition: ERF_SurfaceLayer.H:578
@ MOENG
Moeng 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:581
amrex::Real default_land_surf_moist
Definition: ERF_SurfaceLayer.H:573
amrex::Real surf_temp_flux
Definition: ERF_SurfaceLayer.H:572
amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_SurfaceLayer.H:562
amrex::Real custom_tstar
Definition: ERF_SurfaceLayer.H:577
bool cnk_visc
Definition: ERF_SurfaceLayer.H:580
amrex::Real surf_heating_rate
Definition: ERF_SurfaceLayer.H:571
FluxCalcType flux_type
Definition: ERF_SurfaceLayer.H:553
MoistCalcType moist_type
Definition: ERF_SurfaceLayer.H:555
amrex::Real custom_ustar
Definition: ERF_SurfaceLayer.H:576
amrex::Real default_land_surf_temp
Definition: ERF_SurfaceLayer.H:569
@ SURFACE_TEMPERATURE
Surface temperature specified.
@ HEAT_FLUX
Heat-flux specified.
MOSTAverage m_ma
Definition: ERF_SurfaceLayer.H:587
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
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
175 {
176  // Pointers to the computed averages
177  const auto *const tm_ptr = m_ma.get_average(lev,2); // potential temperature
178  const auto *const qvm_ptr = m_ma.get_average(lev,3); // water vapor mixing ratio
179  const auto *const tvm_ptr = m_ma.get_average(lev,4); // virtual potential temperature
180  const auto *const umm_ptr = m_ma.get_average(lev,5); // horizontal velocity magnitude
181 
182  for (MFIter mfi(*u_star[lev]); mfi.isValid(); ++mfi)
183  {
184  Box gtbx = mfi.growntilebox();
185 
186  auto u_star_arr = u_star[lev]->array(mfi);
187  auto t_star_arr = t_star[lev]->array(mfi);
188  auto q_star_arr = q_star[lev]->array(mfi);
189  auto t_surf_arr = t_surf[lev]->array(mfi);
190  auto q_surf_arr = q_surf[lev]->array(mfi);
191  auto olen_arr = olen[lev]->array(mfi);
192 
193  const auto tm_arr = tm_ptr->array(mfi);
194  const auto tvm_arr = tvm_ptr->array(mfi);
195  const auto qvm_arr = qvm_ptr->array(mfi);
196  const auto umm_arr = umm_ptr->array(mfi);
197  const auto z0_arr = z_0[lev].array();
198 
199  // PBL height if we need to calculate wstar for the Beljaars correction
200  // TODO: can/should we apply this in LES mode?
201  const auto w_star_arr = (m_include_wstar) ? w_star[lev].get()->array(mfi) : Array4<Real> {};
202  const auto pblh_arr = (m_include_wstar) ? pblh[lev].get()->array(mfi) : Array4<Real> {};
203 
204  // Wave properties if they exist
205  const auto Hwave_arr = (m_Hwave_lev[lev]) ? m_Hwave_lev[lev]->array(mfi) : Array4<Real> {};
206  const auto Lwave_arr = (m_Lwave_lev[lev]) ? m_Lwave_lev[lev]->array(mfi) : Array4<Real> {};
207  const auto eta_arr = (m_eddyDiffs_lev[lev]) ? m_eddyDiffs_lev[lev]->array(mfi) : Array4<Real> {};
208 
209  // Land mask array if it exists
210  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
211  Array4<int> {};
212 
213  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
214  {
215  if (( is_land && lmask_arr(i,j,k) == 1) ||
216  (!is_land && lmask_arr(i,j,k) == 0))
217  {
218  most_flux.iterate_flux(i, j, k, max_iters,
219  z0_arr, umm_arr, tm_arr, tvm_arr, qvm_arr,
220  u_star_arr, w_star_arr, // to be updated
221  t_star_arr, q_star_arr, // to be updated
222  t_surf_arr, q_surf_arr, olen_arr, // to be updated
223  pblh_arr, Hwave_arr, Lwave_arr, eta_arr);
224  }
225  });
226  }
227 }
const amrex::MultiFab * get_average(const int &lev, const int &comp) const
Definition: ERF_MOSTAverage.H:101
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
Definition: ERF_SurfaceLayer.H:599
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
Definition: ERF_SurfaceLayer.H:594
amrex::Vector< amrex::FArrayBox > z_0
Definition: ERF_SurfaceLayer.H:582
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
Definition: ERF_SurfaceLayer.H:591
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
Definition: ERF_SurfaceLayer.H:603
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
Definition: ERF_SurfaceLayer.H:589
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
Definition: ERF_SurfaceLayer.H:588
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
Definition: ERF_SurfaceLayer.H:595
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
Definition: ERF_SurfaceLayer.H:602
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
Definition: ERF_SurfaceLayer.H:590
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
Definition: ERF_SurfaceLayer.H:604
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
Definition: ERF_SurfaceLayer.H:592
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
Definition: ERF_SurfaceLayer.H:593

◆ 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 int  RhoQv_comp,
const int  RhoQc_comp,
const int  RhoQr_comp 
)

◆ 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,
int  RhoQv_comp,
int  RhoQc_comp,
int  RhoQr_comp 
)
575 {
576  est.compute_pblh(m_geom[lev],z_phys_cc, pblh[lev].get(),
577  vars[lev][Vars::cons],m_lmask_lev[lev][0],
578  RhoQv_comp, RhoQc_comp, RhoQr_comp);
579 }
@ 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
298 {
299  bool rotate = m_rotate;
300  const int klo = m_geom[lev].Domain().smallEnd(2);
301  const auto& dxInv = m_geom[lev].InvCellSizeArray();
302  for (MFIter mfi(*mfs[0]); mfi.isValid(); ++mfi)
303  {
304  // Get field arrays
305  const auto cons_arr = mfs[Vars::cons]->array(mfi);
306  const auto velx_arr = mfs[Vars::xvel]->array(mfi);
307  const auto vely_arr = mfs[Vars::yvel]->array(mfi);
308  const auto velz_arr = mfs[Vars::zvel]->array(mfi);
309 
310  // Diffusive stress vars
311  auto t13_arr = Tau_lev[TauType::tau13]->array(mfi);
312  auto t31_arr = (Tau_lev[TauType::tau31]) ? Tau_lev[TauType::tau31]->array(mfi) : Array4<Real>{};
313 
314  auto t23_arr = Tau_lev[TauType::tau23]->array(mfi);
315  auto t32_arr = (Tau_lev[TauType::tau32]) ? Tau_lev[TauType::tau32]->array(mfi) : Array4<Real>{};
316 
317  auto hfx3_arr = zheat_flux->array(mfi);
318  auto qfx3_arr = (zqv_flux) ? zqv_flux->array(mfi) : Array4<Real>{};
319 
320  // Rotated stress vars
321  auto t11_arr = (m_rotate) ? Tau_lev[TauType::tau11]->array(mfi) : Array4<Real>{};
322  auto t22_arr = (m_rotate) ? Tau_lev[TauType::tau22]->array(mfi) : Array4<Real>{};
323  auto t33_arr = (m_rotate) ? Tau_lev[TauType::tau33]->array(mfi) : Array4<Real>{};
324  auto t12_arr = (m_rotate) ? Tau_lev[TauType::tau12]->array(mfi) : Array4<Real>{};
325  auto t21_arr = (m_rotate) ? Tau_lev[TauType::tau21]->array(mfi) : Array4<Real>{};
326 
327  auto hfx1_arr = (m_rotate) ? xheat_flux->array(mfi) : Array4<Real>{};
328  auto hfx2_arr = (m_rotate) ? yheat_flux->array(mfi) : Array4<Real>{};
329  auto qfx1_arr = (m_rotate && xqv_flux) ? xqv_flux->array(mfi) : Array4<Real>{};
330  auto qfx2_arr = (m_rotate && yqv_flux) ? yqv_flux->array(mfi) : Array4<Real>{};
331 
332  // Terrain
333  const auto zphys_arr = (z_phys) ? z_phys->const_array(mfi) : Array4<const Real>{};
334 
335  // Get average arrays
336  const auto *const u_mean = m_ma.get_average(lev,0);
337  const auto *const v_mean = m_ma.get_average(lev,1);
338  const auto *const t_mean = m_ma.get_average(lev,2);
339  const auto *const q_mean = m_ma.get_average(lev,3);
340  const auto *const u_mag_mean = m_ma.get_average(lev,5);
341 
342  const auto um_arr = u_mean->array(mfi);
343  const auto vm_arr = v_mean->array(mfi);
344  const auto tm_arr = t_mean->array(mfi);
345  const auto qm_arr = q_mean->array(mfi);
346  const auto umm_arr = u_mag_mean->array(mfi);
347 
348  // Get derived arrays
349  const auto u_star_arr = u_star[lev]->array(mfi);
350  const auto t_star_arr = t_star[lev]->array(mfi);
351  const auto q_star_arr = q_star[lev]->array(mfi);
352  const auto t_surf_arr = t_surf[lev]->array(mfi);
353  const auto q_surf_arr = q_surf[lev]->array(mfi);
354 
355  // Get LSM fluxes
356  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
357  Array4<int> {};
358  auto lsm_flux_arr = (m_lsm_flux_lev[lev][0]) ? m_lsm_flux_lev[lev][0]->array(mfi) :
359  Array4<Real> {};
360 
361  // Rho*Theta flux
362  //============================================================================
363  Box bx = mfi.tilebox();
364  if (bx.smallEnd(2) != klo) { continue; }
365  bx.makeSlab(2,klo);
366  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
367  {
368  Real Tflux = flux_comp.compute_t_flux(i, j, k,
369  cons_arr, velx_arr, vely_arr,
370  umm_arr, tm_arr, u_star_arr,
371  t_star_arr, t_surf_arr);
372 
373  if (rotate) {
374  rotate_scalar_flux(i, j, k, Tflux, dxInv, zphys_arr,
375  hfx1_arr, hfx2_arr, hfx3_arr);
376  } else {
377  hfx3_arr(i,j,klo) = Tflux;
378  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
379  if (is_land && lsm_flux_arr) {
380  lsm_flux_arr(i,j,k) = Tflux;
381  }
382  }
383  });
384 
385  // Rho*Qv flux
386  //============================================================================
387  if (use_moisture) {
388  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
389  {
390  Real Qflux = flux_comp.compute_q_flux(i, j, k,
391  cons_arr, velx_arr, vely_arr,
392  umm_arr, qm_arr, u_star_arr,
393  q_star_arr, q_surf_arr);
394 
395  if (rotate) {
396  rotate_scalar_flux(i, j, k, Qflux, dxInv, zphys_arr,
397  qfx1_arr, qfx2_arr, qfx3_arr);
398  } else {
399  qfx3_arr(i,j,k) = Qflux;
400  }
401  });
402  } // custom
403 
404  // Rho*u flux
405  //============================================================================
406  Box bxx = surroundingNodes(bx,0);
407  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
408  {
409  Real stressx = flux_comp.compute_u_flux(i, j, k,
410  cons_arr, velx_arr, vely_arr,
411  umm_arr, um_arr, u_star_arr);
412 
413  if (rotate) {
414  rotate_stress_tensor(i, j, k, stressx, dxInv, zphys_arr,
415  velx_arr, vely_arr, velz_arr,
416  t11_arr, t22_arr, t33_arr,
417  t12_arr, t21_arr,
418  t13_arr, t31_arr,
419  t23_arr, t32_arr);
420  } else {
421  t13_arr(i,j,k) = stressx;
422  if (t31_arr) { t31_arr(i,j,k) = stressx; }
423  }
424  });
425 
426  // Rho*v flux
427  //============================================================================
428  Box bxy = surroundingNodes(bx,1);
429  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
430  {
431  Real stressy = flux_comp.compute_v_flux(i, j, k,
432  cons_arr, velx_arr, vely_arr,
433  umm_arr, vm_arr, u_star_arr);
434 
435  // NOTE: One stress rotation for ALL the stress components
436  if (!rotate) {
437  t23_arr(i,j,k) = stressy;
438  if (t32_arr) { t32_arr(i,j,k) = stressy; }
439  }
440  });
441  } // mfiter
442 }
@ tau12
Definition: ERF_DataStruct.H:30
@ tau23
Definition: ERF_DataStruct.H:30
@ tau33
Definition: ERF_DataStruct.H:30
@ tau22
Definition: ERF_DataStruct.H:30
@ tau11
Definition: ERF_DataStruct.H:30
@ tau32
Definition: ERF_DataStruct.H:30
@ tau31
Definition: ERF_DataStruct.H:30
@ tau21
Definition: ERF_DataStruct.H:30
@ tau13
Definition: ERF_DataStruct.H:30
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:612
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:633
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
Definition: ERF_SurfaceLayer.H:601
@ 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_tsurf_with_sst_and_tsk()

void SurfaceLayer::fill_tsurf_with_sst_and_tsk ( const int &  lev,
const amrex::Real &  time 
)
447 {
448  int n_times_in_sst = m_sst_lev[lev].size();
449 
450  int n_time_lo, n_time_hi;
451  Real alpha;
452 
453  if (n_times_in_sst > 1) {
454  // Time interpolation
455  Real dT = m_bdy_time_interval;
456  Real time_since_start = time - m_start_bdy_time;
457  int n_time = static_cast<int>( time_since_start / dT);
458  n_time_lo = n_time;
459  n_time_hi = n_time+1;
460  alpha = (time_since_start - n_time * dT) / dT;
461  AMREX_ALWAYS_ASSERT( (n_time >= 0) && (n_time < (m_sst_lev[lev].size()-1)));
462  } else {
463  n_time_lo = 0;
464  n_time_hi = 0;
465  alpha = 1.0;
466  }
467  AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
468 
469  Real oma = 1.0 - alpha;
470 
471  // Define a default land surface temperature if we don't read in tsk
472  Real lst = default_land_surf_temp;
473 
474  bool use_tsk = (m_tsk_lev[lev][0]);
475 
476  // Populate t_surf
477  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
478  {
479  Box gtbx = mfi.growntilebox();
480 
481  auto t_surf_arr = t_surf[lev]->array(mfi);
482 
483  const auto sst_lo_arr = m_sst_lev[lev][n_time_lo]->const_array(mfi);
484  const auto sst_hi_arr = m_sst_lev[lev][n_time_hi]->const_array(mfi);
485 
486  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
487  Array4<int> {};
488 
489  if (use_tsk) {
490  const auto tsk_arr = m_tsk_lev[lev][n_time_lo]->const_array(mfi);
491  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
492  {
493  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
494  if (!is_land) {
495  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
496  + alpha * sst_hi_arr(i,j,k);
497  } else {
498  t_surf_arr(i,j,k) = tsk_arr(i,j,k);
499  }
500  });
501  } else {
502  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
503  {
504  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
505  if (!is_land) {
506  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
507  + alpha * sst_hi_arr(i,j,k);
508  } else {
509  t_surf_arr(i,j,k) = lst;
510  }
511  });
512  }
513  }
514  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
515 }
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
Definition: ERF_SurfaceLayer.H:597
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
Definition: ERF_SurfaceLayer.H:598

◆ get_lmask()

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

◆ get_lsm_tsurf()

void SurfaceLayer::get_lsm_tsurf ( const int &  lev)
519 {
520  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
521  {
522  Box gtbx = mfi.growntilebox();
523 
524  // TODO: LSM does not carry lateral ghost cells.
525  // This copies the valid box into the ghost cells.
526  // Fillboundary is called after this to pick up the
527  // interior ghost and periodic directions. Is there
528  // a better approach?
529  Box vbx = mfi.validbox();
530  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
531  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
532 
533  auto t_surf_arr = t_surf[lev]->array(mfi);
534  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
535  Array4<int> {};
536  const auto lsm_arr = m_lsm_data_lev[lev][0]->const_array(mfi);
537 
538  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
539  {
540  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
541  if (is_land) {
542  int li = amrex::min(amrex::max(i, i_lo), i_hi);
543  int lj = amrex::min(amrex::max(j, j_lo), j_hi);
544  t_surf_arr(i,j,k) = lsm_arr(li,lj,k);
545  }
546  });
547  }
548 }
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
Definition: ERF_SurfaceLayer.H:600

◆ get_mac_avg()

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

◆ get_olen()

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

◆ get_pblh()

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

◆ get_q_star()

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

◆ get_q_surf()

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

◆ get_t_star()

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

◆ get_t_surf()

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

◆ get_u_star()

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

◆ get_w_star()

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

◆ get_z0()

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

◆ get_zref()

amrex::Real SurfaceLayer::get_zref ( )
inline
496 { return m_ma.get_zref(); }
amrex::Real get_zref() const
Definition: ERF_MOSTAverage.H:104
Here is the call graph for this function:

◆ have_variable_sea_roughness()

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

◆ 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
248 {
250  moeng_flux flux_comp;
251  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
252  xheat_flux, yheat_flux, zheat_flux,
253  xqv_flux, yqv_flux, zqv_flux,
254  z_phys, flux_comp);
255  } else if (flux_type == FluxCalcType::DONELAN) {
256  donelan_flux flux_comp;
257  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
258  xheat_flux, yheat_flux, zheat_flux,
259  xqv_flux, yqv_flux, zqv_flux,
260  z_phys, flux_comp);
261  } else if (flux_type == FluxCalcType::ROTATE) {
262  rotate_flux flux_comp;
263  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
264  xheat_flux, yheat_flux, zheat_flux,
265  xqv_flux, yqv_flux, zqv_flux,
266  z_phys, flux_comp);
267  } else {
268  custom_flux flux_comp;
269  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
270  xheat_flux, yheat_flux, zheat_flux,
271  xqv_flux, yqv_flux, zqv_flux,
272  z_phys, flux_comp);
273  }
274 }
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:1869
Definition: ERF_MOSTStress.H:1737
Definition: ERF_MOSTStress.H:1574
Definition: ERF_MOSTStress.H:1983

◆ lmask_min_reduce()

int SurfaceLayer::lmask_min_reduce ( amrex::iMultiFab &  lmask,
const int &  nghost 
)
inline
506  {
507  int lmask_min = amrex::ReduceMin(lmask, nghost, [=] AMREX_GPU_HOST_DEVICE(
508  amrex::Box const& bx, amrex::Array4<int const> const& lm_arr) -> int
509  {
510  int locmin = std::numeric_limits<int>::max();
511  const auto lo = lbound(bx);
512  const auto hi = ubound(bx);
513  for (int j = lo.y; j <= hi.y; ++j) {
514  for (int i = lo.x; i <= hi.x; ++i) {
515  locmin = std::min(locmin, lm_arr(i, j, 0));
516  }
517  }
518  return locmin;
519  });
520 
521  return lmask_min;
522  }

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< amrex::MultiFab * >  lsm_flux,
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
238  {
239  // Update MOST Average
241  Theta_prim, Qv_prim, Qr_prim,
242  z_phys_nd);
243 
244  // Get CC vars
245  amrex::MultiFab& mf = *(mfv[0]);
246 
247  amrex::ParmParse pp("erf");
248 
249  // Do we have a time-varying surface roughness that needs to be saved?
250  if (lev == 0) {
251  const int nghost = 0; // ghost cells not included
252  int lmask_min = lmask_min_reduce(*lmask_lev[0].get(), nghost);
253  amrex::ParallelDescriptor::ReduceIntMin(lmask_min);
254 
255  m_var_z0 = (lmask_min < 1) & (rough_type_sea != RoughCalcType::CONSTANT);
256  if (m_var_z0) {
257  std::string rough_sea_string{"charnock"};
258  pp.query("most.roughness_type_sea", rough_sea_string);
259  amrex::Print() << "Variable sea roughness (type " << rough_sea_string
260  << ")" << std::endl;
261  }
262  }
263 
264  if ((lev == 0) || (lev > nlevs - 1)) {
265  m_Hwave_lev.resize(nlevs);
266  m_Lwave_lev.resize(nlevs);
267  m_eddyDiffs_lev.resize(nlevs);
268 
269  m_lsm_data_lev.resize(nlevs);
270  m_lsm_flux_lev.resize(nlevs);
271 
272  m_sst_lev.resize(nlevs);
273  m_tsk_lev.resize(nlevs);
274  m_lmask_lev.resize(nlevs);
275 
276  // Size the MOST params for all levels
277  z_0.resize(nlevs);
278  u_star.resize(nlevs);
279  w_star.resize(nlevs);
280  t_star.resize(nlevs);
281  q_star.resize(nlevs);
282  t_surf.resize(nlevs);
283  q_surf.resize(nlevs);
284  olen.resize(nlevs);
285  pblh.resize(nlevs);
286  }
287 
288  // Get pointers to SST,TSK and LANDMASK data
289  int nt_tot_sst = sst_lev.size();
290  m_sst_lev[lev].resize(nt_tot_sst);
291  for (int nt(0); nt < nt_tot_sst; ++nt) {
292  m_sst_lev[lev][nt] = sst_lev[nt].get();
293  }
294  int nt_tot_tsk = tsk_lev.size();
295  m_tsk_lev[lev].resize(nt_tot_tsk);
296  for (int nt(0); nt < nt_tot_tsk; ++nt) {
297  m_tsk_lev[lev][nt] = tsk_lev[nt].get();
298  }
299  int nt_tot_lmask = lmask_lev.size();
300  m_lmask_lev[lev].resize(nt_tot_lmask);
301  for (int nt(0); nt < nt_tot_lmask; ++nt) {
302  m_lmask_lev[lev][nt] = lmask_lev[nt].get();
303  }
304 
305  // Get pointers to wave data
306  m_Hwave_lev[lev] = Hwave;
307  m_Lwave_lev[lev] = Lwave;
308  m_eddyDiffs_lev[lev] = eddyDiffs;
309 
310  // Get pointers to LSM data and Fluxes
311  int nvar = lsm_data.size();
312  m_lsm_data_lev[lev].resize(nvar);
313  m_lsm_flux_lev[lev].resize(nvar);
314  for (int n(0); n < nvar; ++n) {
315  m_lsm_data_lev[lev][n] = lsm_data[n];
316  m_lsm_flux_lev[lev][n] = lsm_flux[n];
317  }
318 
319  // Check if there is a user-specified roughness file to be read
320  std::string fname;
321  bool read_z0 = false;
322  if ( (flux_type == FluxCalcType::MOENG) ||
324  read_z0 = pp.query("most.roughness_file_name", fname);
325  }
326 
327  // Attributes for MFs and FABs
328  //--------------------------------------------------------
329  // Create a 2D ba, dm, & ghost cells
330  amrex::BoxArray ba = mf.boxArray();
331  amrex::BoxList bl2d = ba.boxList();
332  for (auto& b : bl2d) {
333  b.setRange(2, 0);
334  }
335  amrex::BoxArray ba2d(std::move(bl2d));
336  const amrex::DistributionMapping& dm = mf.DistributionMap();
337  const int ncomp = 1;
338  amrex::IntVect ng = mf.nGrowVect();
339  ng[2] = 0;
340 
341  // Z0 heights FAB
342  //--------------------------------------------------------
343  amrex::Arena* Arena_Used = amrex::The_Arena();
344 #ifdef AMREX_USE_GPU
345  Arena_Used = amrex::The_Pinned_Arena();
346 #endif
347  amrex::Box bx = amrex::grow(m_geom[lev].Domain(), ng);
348  bx.setSmall(2, 0);
349  bx.setBig(2, 0);
350  z_0[lev].resize(bx, 1, Arena_Used);
351  z_0[lev].setVal<amrex::RunOn::Device>(z0_const);
352  if (read_z0) {
353  read_custom_roughness(lev, fname);
354  }
355 
356  // 2D MFs for U*, T*, T_surf
357  //--------------------------------------------------------
358  u_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
359  u_star[lev]->setVal(1.E34);
360 
361  w_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
362  w_star[lev]->setVal(1.E34);
363 
364  t_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
365  t_star[lev]->setVal(0.0); // default to neutral
366 
367  q_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
368  q_star[lev]->setVal(0.0); // default to dry
369 
370  olen[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
371  olen[lev]->setVal(1.E34);
372 
373  pblh[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
374  pblh[lev]->setVal(1.E34);
375 
376  t_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
377  t_surf[lev]->setVal(default_land_surf_temp);
378 
379  q_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
380  q_surf[lev]->setVal(default_land_surf_moist);
381 
382  // TODO: Do we want lsm_data to have theta at 0 index always?
383  // Do we want an external enum struct for indexing?
384  if (m_sst_lev[lev][0] || m_tsk_lev[lev][0] || m_lsm_data_lev[lev][0]) {
385  // Valid SST, TSK or LSM data; t_surf set before computing fluxes (avoids
386  // extended lambda capture) Note that land temp will be set from m_tsk_lev
387  // while sea temp will be set from m_sst_lev
389  }
390 
391 
392  }
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:70
int lmask_min_reduce(amrex::iMultiFab &lmask, const int &nghost)
Definition: ERF_SurfaceLayer.H:504
void read_custom_roughness(const int &lev, const std::string &fname)
Definition: ERF_SurfaceLayer.cpp:582
@ 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 
)
584 {
585  // Read the file if we are on the coarsest level
586  if (lev==0) {
587  // Only the ioproc reads the file and broadcasts
588  if (ParallelDescriptor::IOProcessor()) {
589  Print()<<"Reading MOST roughness file: "<< fname << std::endl;
590  std::ifstream file(fname);
591  Gpu::HostVector<Real> m_x,m_y,m_z0;
592  Real value1,value2,value3;
593  while(file>>value1>>value2>>value3){
594  m_x.push_back(value1);
595  m_y.push_back(value2);
596  m_z0.push_back(value3);
597  }
598  file.close();
599 
600  // Copy data to the GPU
601  int nnode = m_x.size();
602  Gpu::DeviceVector<Real> d_x(nnode),d_y(nnode),d_z0(nnode);
603  Gpu::copy(Gpu::hostToDevice, m_x.begin(), m_x.end(), d_x.begin());
604  Gpu::copy(Gpu::hostToDevice, m_y.begin(), m_y.end(), d_y.begin());
605  Gpu::copy(Gpu::hostToDevice, m_z0.begin(), m_z0.end(), d_z0.begin());
606  Real* xp = d_x.data();
607  Real* yp = d_y.data();
608  Real* z0p = d_z0.data();
609 
610  // Populate z_phys data
611  Real tol = 1.0e-4;
612  auto dx = m_geom[lev].CellSizeArray();
613  auto ProbLoArr = m_geom[lev].ProbLoArray();
614  int ilo = m_geom[lev].Domain().smallEnd(0);
615  int jlo = m_geom[lev].Domain().smallEnd(1);
616  int klo = 0;
617  int ihi = m_geom[lev].Domain().bigEnd(0);
618  int jhi = m_geom[lev].Domain().bigEnd(1);
619 
620  // Grown box with no z range
621  Box xybx = z_0[lev].box();
622  xybx.setRange(2,0);
623 
624  Array4<Real> const& z0_arr = z_0[lev].array();
625  ParallelFor(xybx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
626  {
627  // Clip indices for ghost-cells
628  int ii = amrex::min(amrex::max(i,ilo),ihi);
629  int jj = amrex::min(amrex::max(j,jlo),jhi);
630 
631  // Location of nodes
632  Real x = ProbLoArr[0] + ii * dx[0];
633  Real y = ProbLoArr[1] + jj * dx[1];
634  int inode = ii + jj * (ihi-ilo+2); // stride is Nx+1
635  if (std::sqrt(std::pow(x-xp[inode],2)+std::pow(y-yp[inode],2)) < tol) {
636  z0_arr(i,j,klo) = z0p[inode];
637  } else {
638  // Unexpected list order, do brute force search
639  Real z0loc = 0.0;
640  bool found = false;
641  for (int n=0; n<nnode; ++n) {
642  Real delta=std::sqrt(std::pow(x-xp[n],2)+std::pow(y-yp[n],2));
643  if (delta < tol) {
644  found = true;
645  z0loc = z0p[n];
646  break;
647  }
648  }
649  AMREX_ASSERT_WITH_MESSAGE(found, "Location read from terrain file does not match the grid!");
650  amrex::ignore_unused(found);
651  z0_arr(i,j,klo) = z0loc;
652  }
653  });
654  } // Is ioproc
655 
656  int ioproc = ParallelDescriptor::IOProcessorNumber();
657  ParallelDescriptor::Barrier();
658  ParallelDescriptor::Bcast(z_0[lev].dataPtr(),z_0[lev].box().numPts(),ioproc);
659  } else {
660  // Create a BC mapper that uses FOEXTRAP at domain bndry
661  Vector<int> bc_lo(3,ERFBCType::foextrap);
662  Vector<int> bc_hi(3,ERFBCType::foextrap);
663  Vector<BCRec> bcr; bcr.push_back(BCRec(bc_lo.data(),bc_hi.data()));
664 
665  // Create ref ratio
666  IntVect ratio;
667  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
668  ratio[idim] = m_geom[lev].Domain().length(idim) / m_geom[0].Domain().length(idim);
669  }
670 
671  // Create interp object and interpolate from the coarsest grid
672  Interpolater* interp = &cell_cons_interp;
673  interp->interp(z_0[0] , 0,
674  z_0[lev], 0,
675  1, z_0[lev].box(),
676  ratio, m_geom[0], m_geom[lev],
677  bcr, 0, 0, RunOn::Gpu);
678  }
679 }
@ m_y
Definition: ERF_DataStruct.H:23
@ m_x
Definition: ERF_DataStruct.H:22
@ foextrap
Definition: ERF_IndexDefines.H:208

Referenced by make_SurfaceLayer_at_level().

Here is the caller graph for this function:

◆ update_fluxes()

void SurfaceLayer::update_fluxes ( const int &  lev,
const amrex::Real &  time,
int  max_iters = 500 
)

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

Parameters
[in]levCurrent level
[in]max_itersmaximum iterations to use
15 {
16  // Update SST data if we have a valid pointer
17  if (m_sst_lev[lev][0]) fill_tsurf_with_sst_and_tsk(lev, time);
18 
19  // TODO: we want 0 index to always be theta?
20  // Update land surface temp if we have a valid pointer
21  if (m_lsm_data_lev[lev][0]) get_lsm_tsurf(lev);
22 
23  // Fill interior ghost cells
24  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
25 
26  // Compute plane averages for all vars (regardless of flux type)
28 
29  // Do we have a constant flux for moisture?
30  bool cons_qflux = ( (moist_type == MoistCalcType::MOISTURE_FLUX) ||
32 
33  // ***************************************************************
34  // Iterate the fluxes if moeng type
35  // First iterate over land -- the only model for surface roughness
36  // over land is RoughCalcType::CONSTANT
37  // ***************************************************************
40  bool is_land = true;
43  surface_flux most_flux(m_ma.get_zref(), surf_temp_flux, surf_moist_flux, cons_qflux);
44  compute_fluxes(lev, max_iters, most_flux, is_land);
45  } else {
46  amrex::Abort("Unknown value for rough_type_land");
47  }
49  update_surf_temp(time);
51  surface_temp most_flux(m_ma.get_zref(), surf_temp_flux, surf_moist_flux, cons_qflux);
52  compute_fluxes(lev, max_iters, most_flux, is_land);
53  } else {
54  amrex::Abort("Unknown value for rough_type_land");
55  }
56  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
60  compute_fluxes(lev, max_iters, most_flux, is_land);
61  } else {
62  amrex::Abort("Unknown value for rough_type_land");
63  }
64  } else {
65  amrex::Abort("Unknown value for theta_type");
66  }
67  } // MOENG -- LAND
68 
69  // ***************************************************************
70  // Iterate the fluxes if moeng type
71  // Next iterate over sea -- the models for surface roughness
72  // over sea are CHARNOCK, DONELAN, MODIFIED_CHARNOCK or WAVE_COUPLED
73  // ***************************************************************
76  bool is_land = false;
81  cnk_a, cnk_visc, cons_qflux);
82  compute_fluxes(lev, max_iters, most_flux, is_land);
86  depth, cons_qflux);
87  compute_fluxes(lev, max_iters, most_flux, is_land);
88  } else if (rough_type_sea == RoughCalcType::DONELAN) {
91  cons_qflux);
92  compute_fluxes(lev, max_iters, most_flux, is_land);
96  cons_qflux);
97  compute_fluxes(lev, max_iters, most_flux, is_land);
98  } else {
99  amrex::Abort("Unknown value for rough_type_sea");
100  }
101 
103  update_surf_temp(time);
105  surface_temp_charnock most_flux(m_ma.get_zref(),
107  cnk_a, cnk_visc, cons_qflux);
108  compute_fluxes(lev, max_iters, most_flux, is_land);
112  depth, cons_qflux);
113  compute_fluxes(lev, max_iters, most_flux, is_land);
114  } else if (rough_type_sea == RoughCalcType::DONELAN) {
115  surface_temp_donelan most_flux(m_ma.get_zref(),
117  cons_qflux);
118  compute_fluxes(lev, max_iters, most_flux, is_land);
122  cons_qflux);
123  compute_fluxes(lev, max_iters, most_flux, is_land);
124  } else {
125  amrex::Abort("Unknown value for rough_type_sea");
126  }
127 
128  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
131  adiabatic_charnock most_flux(m_ma.get_zref(),
133  cnk_a, cnk_visc);
134  compute_fluxes(lev, max_iters, most_flux, is_land);
138  depth);
139  compute_fluxes(lev, max_iters, most_flux, is_land);
140  } else if (rough_type_sea == RoughCalcType::DONELAN) {
142  compute_fluxes(lev, max_iters, most_flux, is_land);
145  compute_fluxes(lev, max_iters, most_flux, is_land);
146  } else {
147  amrex::Abort("Unknown value for rough_type_sea");
148  }
149  } else {
150  amrex::Abort("Unknown value for theta_type");
151  }
152 
153  } // MOENG -- SEA
154 
156  u_star[lev]->setVal(custom_ustar);
157  t_star[lev]->setVal(custom_tstar);
158  q_star[lev]->setVal(custom_qstar);
159  }
160 }
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:679
void update_surf_temp(const amrex::Real &time)
Definition: ERF_SurfaceLayer.H:454
void get_lsm_tsurf(const int &lev)
Definition: ERF_SurfaceLayer.cpp:518
void compute_fluxes(const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
Definition: ERF_SurfaceLayer.cpp:171
void fill_tsurf_with_sst_and_tsk(const int &lev, const amrex::Real &time)
Definition: ERF_SurfaceLayer.cpp:445
Definition: ERF_MOSTStress.H:141
Definition: ERF_MOSTStress.H:283
Definition: ERF_MOSTStress.H:216
Definition: ERF_MOSTStress.H:350
Definition: ERF_MOSTStress.H:91
Definition: ERF_MOSTStress.H:526
Definition: ERF_MOSTStress.H:752
Definition: ERF_MOSTStress.H:643
Definition: ERF_MOSTStress.H:858
Definition: ERF_MOSTStress.H:423
Definition: ERF_MOSTStress.H:1087
Definition: ERF_MOSTStress.H:1333
Definition: ERF_MOSTStress.H:1214
Definition: ERF_MOSTStress.H:1449
Definition: ERF_MOSTStress.H:973

◆ 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
471  {
472  m_ma.update_field_ptrs(lev, vars_old, Theta_prim, Qv_prim, Qr_prim);
473  }
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:246
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 int  RhoQv_comp,
const int  RhoQc_comp,
const int  RhoQr_comp 
)
557 {
559  MYNNPBLH estimator;
560  compute_pblh(lev, vars, z_phys_cc, estimator, RhoQv_comp, RhoQc_comp, RhoQr_comp);
562  amrex::Error("YSU/MRF PBLH calc not implemented yet");
563  }
564 }
void compute_pblh(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const int RhoQv_comp, const int RhoQc_comp, const int RhoQr_comp)
Definition: ERF_PBLHeight.H:8

◆ update_surf_temp()

void SurfaceLayer::update_surf_temp ( const amrex::Real &  time)
inline
455  {
456  if (surf_heating_rate != 0) {
457  int nlevs = m_geom.size();
458  for (int lev = 0; lev < nlevs; lev++) {
459  t_surf[lev]->setVal(surf_temp + surf_heating_rate * time);
460  amrex::Print() << "Surface temp at t=" << time << ": "
461  << surf_temp + surf_heating_rate * time << std::endl;
462  }
463  }
464  }

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_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}

◆ m_bdy_time_interval

amrex::Real SurfaceLayer::m_bdy_time_interval
private

◆ m_eddyDiffs_lev

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

◆ m_geom

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

◆ m_Hwave_lev

amrex::Vector<amrex::MultiFab*> SurfaceLayer::m_Hwave_lev
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_flux_lev

amrex::Vector<amrex::Vector<amrex::MultiFab*> > SurfaceLayer::m_lsm_flux_lev
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_bdy_time

amrex::Real SurfaceLayer::m_start_bdy_time
private

◆ m_tsk_lev

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

◆ 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

◆ rough_type_land

RoughCalcType SurfaceLayer::rough_type_land {RoughCalcType::CONSTANT}

Referenced by SurfaceLayer().

◆ rough_type_sea

RoughCalcType SurfaceLayer::rough_type_sea {RoughCalcType::CHARNOCK}

◆ 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_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::FArrayBox> SurfaceLayer::z_0
private

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