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 , 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 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, const amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, 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 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)
 
amrex::MultiFab * get_q_surf (const int &lev)
 
amrex::Real get_zref ()
 
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)
 
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_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}
 
bool use_moisture
 
bool m_has_lsm_tsurf = false
 
int m_lsm_tsurf_indx = -1
 
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
 

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.

531  {
532  MOENG = 0, ///< Moeng functional form
533  DONELAN, ///< Donelan functional form
534  CUSTOM, ///< Custom constant flux functional form
535  ROTATE ///< Terrain rotation flux functional form
536  };

◆ MoistCalcType

Enumerator
ADIABATIC 
MOISTURE_FLUX 

Qv-flux specified.

SURFACE_MOISTURE 

Surface Qv specified.

544  {
545  ADIABATIC = 0,
546  MOISTURE_FLUX, ///< Qv-flux specified
547  SURFACE_MOISTURE ///< Surface Qv specified
548  };

◆ PBLHeightCalcType

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

◆ RoughCalcType

Enumerator
CONSTANT 

Constant z0.

CHARNOCK 
MODIFIED_CHARNOCK 
DONELAN 
WAVE_COUPLED 
550  {
551  CONSTANT = 0, ///< Constant z0
552  CHARNOCK,
553  MODIFIED_CHARNOCK,
554  DONELAN,
555  WAVE_COUPLED
556  };

◆ ThetaCalcType

Enumerator
ADIABATIC 
HEAT_FLUX 

Heat-flux specified.

SURFACE_TEMPERATURE 

Surface temperature specified.

538  {
539  ADIABATIC = 0,
540  HEAT_FLUX, ///< Heat-flux specified
541  SURFACE_TEMPERATURE ///< Surface temperature specified
542  };

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  bdy_time_interval = 0.0 
)
inlineexplicit
42  : m_geom(geom),
43  m_rotate(use_rot_surface_flux),
44  m_bdy_time_interval(bdy_time_interval),
45  m_ma(geom, (z_phys_nd[0] != nullptr), a_pp_prefix, a_terrain_type)
46  {
47  // We have a moisture model if Qv_prim is a valid pointer
48  use_moisture = (Qv_prim[0].get());
49 
50  // Get roughness
51  amrex::ParmParse pp("erf");
52  pp.query("most.z0", z0_const);
53 
54  // Specify how to compute the flux
55  if (use_rot_surface_flux) {
57  } else {
58  std::string flux_string_in;
59  std::string flux_string{"moeng"};
60  auto read_flux = pp.query("surface_layer.flux_type", flux_string_in);
61  if (read_flux) {
62  flux_string = amrex::toLower(flux_string_in);
63  }
64  if (flux_string == "donelan") {
66  } else if (flux_string == "moeng") {
68  } else if (flux_string == "custom") {
70  } else {
71  amrex::Abort("Undefined MOST flux type!");
72  }
73  }
74 
75  // Include w* to handle free convection (Beljaars 1995, QJRMS)
76  pp.query("most.include_wstar", m_include_wstar);
77 
78  std::string pblh_string_in;
79  std::string pblh_string{"none"};
80  auto read_pblh = pp.query("most.pblh_calc", pblh_string_in);
81  if (read_pblh) {
82  pblh_string = amrex::toLower(pblh_string_in);
83  }
84  if (pblh_string == "none") {
86  } else if (pblh_string == "mynn25") {
88  } else if (pblh_string == "mynnedmf") {
90  } else if (pblh_string == "ysu") {
92  } else if (pblh_string == "mrf") {
94  } else {
95  amrex::Abort("Undefined PBLH calc type!");
96  }
97 
98  // Get surface temperature
99  auto erf_st = pp.query("most.surf_temp", surf_temp);
100  if (erf_st) { default_land_surf_temp = surf_temp; }
101 
102  // Get surface moisture
103  bool erf_sq = false;
104  if (use_moisture) { erf_sq = pp.query("most.surf_moist", surf_moist); }
105  if (erf_sq) { default_land_surf_moist = surf_moist; }
106 
107  // Custom type user must specify the fluxes
111  pp.query("most.ustar", custom_ustar);
112  pp.query("most.tstar", custom_tstar);
113  pp.query("most.qstar", custom_qstar);
114  pp.query("most.rhosurf", custom_rhosurf);
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:230
ThetaCalcType theta_type
Definition: ERF_SurfaceLayer.H:561
bool m_include_wstar
Definition: ERF_SurfaceLayer.H:573
bool m_rotate
Definition: ERF_SurfaceLayer.H:570
PBLHeightCalcType pblh_type
Definition: ERF_SurfaceLayer.H:565
bool use_moisture
Definition: ERF_SurfaceLayer.H:593
RoughCalcType rough_type_land
Definition: ERF_SurfaceLayer.H:563
amrex::Real z0_const
Definition: ERF_SurfaceLayer.H:574
amrex::Real cnk_a
Definition: ERF_SurfaceLayer.H:587
amrex::Real surf_temp
Definition: ERF_SurfaceLayer.H:576
amrex::Real surf_moist_flux
Definition: ERF_SurfaceLayer.H:581
RoughCalcType rough_type_sea
Definition: ERF_SurfaceLayer.H:564
amrex::Real surf_moist
Definition: ERF_SurfaceLayer.H:580
amrex::Real m_bdy_time_interval
Definition: ERF_SurfaceLayer.H:571
amrex::Real custom_qstar
Definition: ERF_SurfaceLayer.H:584
amrex::Real custom_rhosurf
Definition: ERF_SurfaceLayer.H:585
@ 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:589
amrex::Real default_land_surf_moist
Definition: ERF_SurfaceLayer.H:579
amrex::Real surf_temp_flux
Definition: ERF_SurfaceLayer.H:578
amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_SurfaceLayer.H:569
amrex::Real custom_tstar
Definition: ERF_SurfaceLayer.H:583
bool cnk_visc
Definition: ERF_SurfaceLayer.H:588
amrex::Real surf_heating_rate
Definition: ERF_SurfaceLayer.H:577
FluxCalcType flux_type
Definition: ERF_SurfaceLayer.H:560
MoistCalcType moist_type
Definition: ERF_SurfaceLayer.H:562
amrex::Real custom_ustar
Definition: ERF_SurfaceLayer.H:582
amrex::Real default_land_surf_temp
Definition: ERF_SurfaceLayer.H:575
@ SURFACE_TEMPERATURE
Surface temperature specified.
@ HEAT_FLUX
Heat-flux specified.
MOSTAverage m_ma
Definition: ERF_SurfaceLayer.H:597
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
193 {
194  // Pointers to the computed averages
195  const auto *const tm_ptr = m_ma.get_average(lev,2); // potential temperature
196  const auto *const qvm_ptr = m_ma.get_average(lev,3); // water vapor mixing ratio
197  const auto *const tvm_ptr = m_ma.get_average(lev,4); // virtual potential temperature
198  const auto *const umm_ptr = m_ma.get_average(lev,5); // horizontal velocity magnitude
199 
200  for (MFIter mfi(*u_star[lev]); mfi.isValid(); ++mfi)
201  {
202  Box gtbx = mfi.growntilebox();
203 
204  auto u_star_arr = u_star[lev]->array(mfi);
205  auto t_star_arr = t_star[lev]->array(mfi);
206  auto q_star_arr = q_star[lev]->array(mfi);
207  auto t_surf_arr = t_surf[lev]->array(mfi);
208  auto q_surf_arr = q_surf[lev]->array(mfi);
209  auto olen_arr = olen[lev]->array(mfi);
210 
211  const auto tm_arr = tm_ptr->array(mfi);
212  const auto tvm_arr = tvm_ptr->array(mfi);
213  const auto qvm_arr = qvm_ptr->array(mfi);
214  const auto umm_arr = umm_ptr->array(mfi);
215  const auto z0_arr = z_0[lev].array(mfi);
216 
217  // PBL height if we need to calculate wstar for the Beljaars correction
218  // TODO: can/should we apply this in LES mode?
219  const auto w_star_arr = (m_include_wstar) ? w_star[lev].get()->array(mfi) : Array4<Real> {};
220  const auto pblh_arr = (m_include_wstar) ? pblh[lev].get()->array(mfi) : Array4<Real> {};
221 
222  // Wave properties if they exist
223  const auto Hwave_arr = (m_Hwave_lev[lev]) ? m_Hwave_lev[lev]->array(mfi) : Array4<Real> {};
224  const auto Lwave_arr = (m_Lwave_lev[lev]) ? m_Lwave_lev[lev]->array(mfi) : Array4<Real> {};
225  const auto eta_arr = (m_eddyDiffs_lev[lev]) ? m_eddyDiffs_lev[lev]->array(mfi) : Array4<Real> {};
226 
227  // Land mask array if it exists
228  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
229  Array4<int> {};
230 
231  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
232  {
233  if (( is_land && lmask_arr(i,j,k) == 1) ||
234  (!is_land && lmask_arr(i,j,k) == 0))
235  {
236  most_flux.iterate_flux(i, j, k, max_iters,
237  z0_arr, umm_arr, tm_arr, tvm_arr, qvm_arr,
238  u_star_arr, w_star_arr, // to be updated
239  t_star_arr, q_star_arr, // to be updated
240  t_surf_arr, q_surf_arr, olen_arr, // to be updated
241  pblh_arr, Hwave_arr, Lwave_arr, eta_arr);
242  }
243  });
244  }
245 }
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:609
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
Definition: ERF_SurfaceLayer.H:604
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
Definition: ERF_SurfaceLayer.H:601
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
Definition: ERF_SurfaceLayer.H:615
amrex::Vector< amrex::MultiFab > z_0
Definition: ERF_SurfaceLayer.H:590
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
Definition: ERF_SurfaceLayer.H:599
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
Definition: ERF_SurfaceLayer.H:598
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
Definition: ERF_SurfaceLayer.H:605
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
Definition: ERF_SurfaceLayer.H:614
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
Definition: ERF_SurfaceLayer.H:600
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
Definition: ERF_SurfaceLayer.H:616
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
Definition: ERF_SurfaceLayer.H:602
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
Definition: ERF_SurfaceLayer.H:603

◆ 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 
)
669 {
670  est.compute_pblh(m_geom[lev],z_phys_cc, pblh[lev].get(),
671  vars[lev][Vars::cons],m_lmask_lev[lev][0],
672  moisture_indices);
673 }
@ 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
316 {
317  bool rotate = m_rotate;
318  const int klo = m_geom[lev].Domain().smallEnd(2);
319  const auto& dxInv = m_geom[lev].InvCellSizeArray();
320  for (MFIter mfi(*mfs[0]); mfi.isValid(); ++mfi)
321  {
322  // Get field arrays
323  const auto cons_arr = mfs[Vars::cons]->array(mfi);
324  const auto velx_arr = mfs[Vars::xvel]->array(mfi);
325  const auto vely_arr = mfs[Vars::yvel]->array(mfi);
326  const auto velz_arr = mfs[Vars::zvel]->array(mfi);
327 
328  // Diffusive stress vars
329  auto t13_arr = Tau_lev[TauType::tau13]->array(mfi);
330  auto t31_arr = (Tau_lev[TauType::tau31]) ? Tau_lev[TauType::tau31]->array(mfi) : Array4<Real>{};
331 
332  auto t23_arr = Tau_lev[TauType::tau23]->array(mfi);
333  auto t32_arr = (Tau_lev[TauType::tau32]) ? Tau_lev[TauType::tau32]->array(mfi) : Array4<Real>{};
334 
335  auto hfx3_arr = zheat_flux->array(mfi);
336  auto qfx3_arr = (zqv_flux) ? zqv_flux->array(mfi) : Array4<Real>{};
337 
338  // Rotated stress vars
339  auto t11_arr = (m_rotate) ? Tau_lev[TauType::tau11]->array(mfi) : Array4<Real>{};
340  auto t22_arr = (m_rotate) ? Tau_lev[TauType::tau22]->array(mfi) : Array4<Real>{};
341  auto t33_arr = (m_rotate) ? Tau_lev[TauType::tau33]->array(mfi) : Array4<Real>{};
342  auto t12_arr = (m_rotate) ? Tau_lev[TauType::tau12]->array(mfi) : Array4<Real>{};
343  auto t21_arr = (m_rotate) ? Tau_lev[TauType::tau21]->array(mfi) : Array4<Real>{};
344 
345  auto hfx1_arr = (m_rotate) ? xheat_flux->array(mfi) : Array4<Real>{};
346  auto hfx2_arr = (m_rotate) ? yheat_flux->array(mfi) : Array4<Real>{};
347  auto qfx1_arr = (m_rotate && xqv_flux) ? xqv_flux->array(mfi) : Array4<Real>{};
348  auto qfx2_arr = (m_rotate && yqv_flux) ? yqv_flux->array(mfi) : Array4<Real>{};
349 
350  // Terrain
351  const auto zphys_arr = (z_phys) ? z_phys->const_array(mfi) : Array4<const Real>{};
352 
353  // Get average arrays
354  const auto *const u_mean = m_ma.get_average(lev,0);
355  const auto *const v_mean = m_ma.get_average(lev,1);
356  const auto *const t_mean = m_ma.get_average(lev,2);
357  const auto *const q_mean = m_ma.get_average(lev,3);
358  const auto *const u_mag_mean = m_ma.get_average(lev,5);
359 
360  const auto um_arr = u_mean->array(mfi);
361  const auto vm_arr = v_mean->array(mfi);
362  const auto tm_arr = t_mean->array(mfi);
363  const auto qm_arr = q_mean->array(mfi);
364  const auto umm_arr = u_mag_mean->array(mfi);
365 
366  // Get derived arrays
367  const auto u_star_arr = u_star[lev]->array(mfi);
368  const auto t_star_arr = t_star[lev]->array(mfi);
369  const auto q_star_arr = q_star[lev]->array(mfi);
370  const auto t_surf_arr = t_surf[lev]->array(mfi);
371  const auto q_surf_arr = q_surf[lev]->array(mfi);
372 
373  // Get LSM fluxes
374  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
375  Array4<int> {};
376  auto lsm_t_flux_arr = Array4<Real> {};
377  auto lsm_q_flux_arr = Array4<Real> {};
378  auto lsm_tau13_arr = Array4<Real> {};
379  auto lsm_tau23_arr = Array4<Real> {};
380  for (int n(0); n<m_lsm_flux_lev[lev].size(); ++n) {
381  if (toLower(m_lsm_flux_name[n]) == "t_flux") { lsm_t_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
382  if (toLower(m_lsm_flux_name[n]) == "q_flux") { lsm_q_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
383  if (toLower(m_lsm_flux_name[n]) == "tau13") { lsm_tau13_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
384  if (toLower(m_lsm_flux_name[n]) == "tau23") { lsm_tau23_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
385  }
386 
387 
388  // Rho*Theta flux
389  //============================================================================
390  Box bx = mfi.tilebox();
391  if (bx.smallEnd(2) != klo) { continue; }
392  bx.makeSlab(2,klo);
393  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
394  {
395  // Valid theta flux from LSM and over land
396  Real Tflux;
397  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
398  if (lsm_t_flux_arr && is_land) {
399  Tflux = lsm_t_flux_arr(i,j,k);
400  } else {
401  Tflux = flux_comp.compute_t_flux(i, j, k,
402  cons_arr, velx_arr, vely_arr,
403  umm_arr, tm_arr, u_star_arr,
404  t_star_arr, t_surf_arr);
405 
406  }
407 
408  // Do scalar flux rotations?
409  if (rotate) {
410  rotate_scalar_flux(i, j, k, Tflux, dxInv, zphys_arr,
411  hfx1_arr, hfx2_arr, hfx3_arr);
412  } else {
413  hfx3_arr(i,j,klo) = Tflux;
414  }
415  });
416 
417  // Rho*Qv flux
418  //============================================================================
419  if (use_moisture) {
420  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
421  {
422  // Valid qv flux from LSM and over land
423  Real Qflux;
424  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
425  if (lsm_q_flux_arr && is_land) {
426  Qflux = lsm_q_flux_arr(i,j,k);
427  } else {
428  Qflux = flux_comp.compute_q_flux(i, j, k,
429  cons_arr, velx_arr, vely_arr,
430  umm_arr, qm_arr, u_star_arr,
431  q_star_arr, q_surf_arr);
432  }
433 
434  // Do scalar flux rotations?
435  if (rotate) {
436  rotate_scalar_flux(i, j, k, Qflux, dxInv, zphys_arr,
437  qfx1_arr, qfx2_arr, qfx3_arr);
438  } else {
439  qfx3_arr(i,j,k) = Qflux;
440  }
441  });
442  } // custom
443 
444  if (!rotate) {
445  // Rho*u flux
446  //============================================================================
447  Box bxx = surroundingNodes(bx,0);
448  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
449  {
450  // Valid tau13 from LSM and over land
451  Real stressx;
452  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
453  if (lsm_tau13_arr && is_land) {
454  stressx = lsm_tau13_arr(i,j,k);
455  } else {
456  stressx = flux_comp.compute_u_flux(i, j, k,
457  cons_arr, velx_arr, vely_arr,
458  umm_arr, um_arr, u_star_arr);
459  }
460 
461  t13_arr(i,j,k) = stressx;
462  if (t31_arr) { t31_arr(i,j,k) = stressx; }
463  });
464 
465  // Rho*v flux
466  //============================================================================
467  Box bxy = surroundingNodes(bx,1);
468  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
469  {
470  // Valid tau13 from LSM and over land
471  Real stressy;
472  int is_land = (lmask_arr) ? lmask_arr(i,j,klo) : 1;
473  if (lsm_tau23_arr && is_land) {
474  stressy = lsm_tau23_arr(i,j,k);
475  } else {
476  stressy = flux_comp.compute_v_flux(i, j, k,
477  cons_arr, velx_arr, vely_arr,
478  umm_arr, vm_arr, u_star_arr);
479  }
480 
481  t23_arr(i,j,k) = stressy;
482  if (t32_arr) { t32_arr(i,j,k) = stressy; }
483  });
484  } else {
485  // All fluxes with rotation
486  //============================================================================
487  Box bxxy = convert(bx, IntVect(1,1,0));
488  ParallelFor(bxxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
489  {
490  Real stresst = flux_comp.compute_u_flux(i, j, k,
491  cons_arr, velx_arr, vely_arr,
492  umm_arr, um_arr, u_star_arr);
493  rotate_stress_tensor(i, j, k, stresst, dxInv, zphys_arr,
494  velx_arr, vely_arr, velz_arr,
495  t11_arr, t22_arr, t33_arr,
496  t12_arr, t21_arr,
497  t13_arr, t31_arr,
498  t23_arr, t32_arr);
499  });
500  }
501  } // mfiter
502 }
@ 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::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: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:611
amrex::Vector< std::string > m_lsm_flux_name
Definition: ERF_SurfaceLayer.H:613
@ 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 
)
580 {
581  // NOTE: We have already tested a moisture model exists
582 
583  // Populate q_surf with qsat over water
584  auto dz = m_geom[lev].CellSize(2);
585  for (MFIter mfi(*q_surf[lev]); mfi.isValid(); ++mfi)
586  {
587  Box gtbx = mfi.growntilebox();
588 
589  auto t_surf_arr = t_surf[lev]->array(mfi);
590  auto q_surf_arr = q_surf[lev]->array(mfi);
591  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
592  Array4<int> {};
593  const auto cons_arr = cons_in.const_array(mfi);
594  const auto z_arr = (z_phys_nd) ? z_phys_nd->const_array(mfi) :
595  Array4<const Real> {};
596 
597  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
598  {
599  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
600  if (!is_land) {
601  auto deltaZ = (z_arr) ? Compute_Zrel_AtCellCenter(i,j,k,z_arr) :
602  0.5*dz;
603  auto Rho = cons_arr(i,j,k,Rho_comp);
604  auto RTh = cons_arr(i,j,k,RhoTheta_comp);
605  auto Qv = cons_arr(i,j,k,RhoQ1_comp) / Rho;
606  auto P_cc = getPgivenRTh(RTh, Qv);
607  P_cc += Rho*CONST_GRAV*deltaZ;
608  P_cc *= 0.01;
609  erf_qsatw(t_surf_arr(i,j,k), P_cc, q_surf_arr(i,j,k));
610  }
611  });
612  }
613  q_surf[lev]->FillBoundary(m_geom[lev].periodicity());
614 }
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:163
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:390
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 
)
507 {
508  int n_times_in_sst = m_sst_lev[lev].size();
509 
510  int n_time_lo, n_time_hi;
511  Real alpha;
512 
513  if (n_times_in_sst > 1) {
514  // Time interpolation
516  int n_time = static_cast<int>( time / dT);
517  n_time_lo = n_time;
518  n_time_hi = n_time+1;
519  alpha = (time - n_time * dT) / dT;
520  AMREX_ALWAYS_ASSERT( (n_time >= 0) && (n_time < (m_sst_lev[lev].size()-1)));
521  } else {
522  n_time_lo = 0;
523  n_time_hi = 0;
524  alpha = 1.0;
525  }
526  AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
527 
528  Real oma = 1.0 - alpha;
529 
530  // Define a default land surface temperature if we don't read in tsk
532 
533  bool use_tsk = (m_tsk_lev[lev][0]);
534 
535  // Populate t_surf
536  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
537  {
538  Box gtbx = mfi.growntilebox();
539 
540  auto t_surf_arr = t_surf[lev]->array(mfi);
541 
542  const auto sst_lo_arr = m_sst_lev[lev][n_time_lo]->const_array(mfi);
543  const auto sst_hi_arr = m_sst_lev[lev][n_time_hi]->const_array(mfi);
544 
545  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
546  Array4<int> {};
547 
548  if (use_tsk) {
549  const auto tsk_arr = m_tsk_lev[lev][n_time_lo]->const_array(mfi);
550  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
551  {
552  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
553  if (!is_land) {
554  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
555  + alpha * sst_hi_arr(i,j,k);
556  } else {
557  t_surf_arr(i,j,k) = tsk_arr(i,j,k);
558  }
559  });
560  } else {
561  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
562  {
563  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
564  if (!is_land) {
565  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
566  + alpha * sst_hi_arr(i,j,k);
567  } else {
568  t_surf_arr(i,j,k) = lst;
569  }
570  });
571  }
572  }
573  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
574 }
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
Definition: ERF_SurfaceLayer.H:607
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
Definition: ERF_SurfaceLayer.H:608

◆ get_lmask()

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

◆ get_lsm_tsurf()

void SurfaceLayer::get_lsm_tsurf ( const int &  lev)
618 {
619  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
620  {
621  Box gtbx = mfi.growntilebox();
622 
623  // NOTE: LSM does not carry lateral ghost cells.
624  // This copies the valid box into the ghost cells.
625  // Fillboundary is called after this to pick up the
626  // interior ghost and periodic directions.
627  Box vbx = mfi.validbox();
628  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
629  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
630 
631  auto t_surf_arr = t_surf[lev]->array(mfi);
632  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
633  Array4<int> {};
634  const auto lsm_arr = m_lsm_data_lev[lev][m_lsm_tsurf_indx]->const_array(mfi);
635 
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) {
640  int li = amrex::min(amrex::max(i, i_lo), i_hi);
641  int lj = amrex::min(amrex::max(j, j_lo), j_hi);
642  t_surf_arr(i,j,k) = lsm_arr(li,lj,k);
643  }
644  });
645  }
646 }
int m_lsm_tsurf_indx
Definition: ERF_SurfaceLayer.H:595
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
Definition: ERF_SurfaceLayer.H:610

◆ get_mac_avg()

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

◆ get_olen()

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

◆ get_pblh()

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

◆ get_q_star()

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

◆ get_q_surf()

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

◆ get_t_star()

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

◆ get_t_surf()

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

◆ get_u_star()

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

◆ get_w_star()

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

◆ get_z0()

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

◆ get_zref()

amrex::Real SurfaceLayer::get_zref ( )
inline
503 { 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
507 { return m_var_z0; }
bool m_var_z0
Definition: ERF_SurfaceLayer.H:591

◆ 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
266 {
268  moeng_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  } else if (flux_type == FluxCalcType::DONELAN) {
274  donelan_flux flux_comp;
275  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
276  xheat_flux, yheat_flux, zheat_flux,
277  xqv_flux, yqv_flux, zqv_flux,
278  z_phys, flux_comp);
279  } else if (flux_type == FluxCalcType::ROTATE) {
280  rotate_flux flux_comp;
281  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
282  xheat_flux, yheat_flux, zheat_flux,
283  xqv_flux, yqv_flux, zqv_flux,
284  z_phys, flux_comp);
285  } else {
286  custom_flux flux_comp(specified_rho_surf);
287  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
288  xheat_flux, yheat_flux, zheat_flux,
289  xqv_flux, yqv_flux, zqv_flux,
290  z_phys, flux_comp);
291  }
292 }
bool specified_rho_surf
Definition: ERF_SurfaceLayer.H:586
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:2061
Definition: ERF_MOSTStress.H:1929
Definition: ERF_MOSTStress.H:1766
Definition: ERF_MOSTStress.H:2178

◆ lmask_min_reduce()

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

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
240  {
241  // Update MOST Average
243  Theta_prim, Qv_prim, Qr_prim,
244  z_phys_nd);
245 
246  // Get CC vars
247  amrex::MultiFab& mf = *(mfv[0]);
248 
249  amrex::ParmParse pp("erf");
250 
251  // Do we have a time-varying surface roughness that needs to be saved?
252  if (lev == 0) {
253  const int nghost = 0; // ghost cells not included
254  int lmask_min = lmask_min_reduce(*lmask_lev[0].get(), nghost);
255  amrex::ParallelDescriptor::ReduceIntMin(lmask_min);
256 
257  m_var_z0 = (lmask_min < 1) & (rough_type_sea != RoughCalcType::CONSTANT);
258  if (m_var_z0) {
259  std::string rough_sea_string{"charnock"};
260  pp.query("most.roughness_type_sea", rough_sea_string);
261  amrex::Print() << "Variable sea roughness (type " << rough_sea_string
262  << ")" << std::endl;
263  }
264  }
265 
266  if ((lev == 0) || (lev > nlevs - 1)) {
267  m_Hwave_lev.resize(nlevs);
268  m_Lwave_lev.resize(nlevs);
269  m_eddyDiffs_lev.resize(nlevs);
270 
271  m_lsm_data_lev.resize(nlevs);
272  m_lsm_flux_lev.resize(nlevs);
273 
274  m_sst_lev.resize(nlevs);
275  m_tsk_lev.resize(nlevs);
276  m_lmask_lev.resize(nlevs);
277 
278  // Size the MOST params for all levels
279  z_0.resize(nlevs);
280  u_star.resize(nlevs);
281  w_star.resize(nlevs);
282  t_star.resize(nlevs);
283  q_star.resize(nlevs);
284  t_surf.resize(nlevs);
285  q_surf.resize(nlevs);
286  olen.resize(nlevs);
287  pblh.resize(nlevs);
288  }
289 
290  // Get pointers to SST,TSK and LANDMASK data
291  int nt_tot_sst = sst_lev.size();
292  m_sst_lev[lev].resize(nt_tot_sst);
293  for (int nt(0); nt < nt_tot_sst; ++nt) {
294  m_sst_lev[lev][nt] = sst_lev[nt].get();
295  }
296  int nt_tot_tsk = static_cast<int>(tsk_lev.size());
297  m_tsk_lev[lev].resize(nt_tot_tsk);
298  for (int nt(0); nt < nt_tot_tsk; ++nt) {
299  m_tsk_lev[lev][nt] = tsk_lev[nt].get();
300  }
301  int nt_tot_lmask = static_cast<int>(lmask_lev.size());
302  m_lmask_lev[lev].resize(nt_tot_lmask);
303  for (int nt(0); nt < nt_tot_lmask; ++nt) {
304  m_lmask_lev[lev][nt] = lmask_lev[nt].get();
305  }
306 
307  // Get pointers to wave data
308  m_Hwave_lev[lev] = Hwave;
309  m_Lwave_lev[lev] = Lwave;
310  m_eddyDiffs_lev[lev] = eddyDiffs;
311 
312  // Get pointers to LSM data and Fluxes
313  int ndata = static_cast<int>(lsm_data.size());
314  int nflux = static_cast<int>(lsm_flux.size());
315  m_lsm_data_name.resize(ndata);
316  m_lsm_data_lev[lev].resize(ndata);
317  m_lsm_flux_name.resize(nflux);
318  m_lsm_flux_lev[lev].resize(nflux);
319  for (int n(0); n < ndata; ++n) {
320  m_lsm_data_name[n] = lsm_data_name[n];
321  m_lsm_data_lev[lev][n] = lsm_data[n];
322  if (amrex::toLower(lsm_data_name[n]) == "theta") {
323  m_has_lsm_tsurf = true;
324  m_lsm_tsurf_indx = n;
325  }
326  }
327  for (int n(0); n < nflux; ++n) {
328  m_lsm_flux_name[n] = lsm_flux_name[n];
329  m_lsm_flux_lev[lev][n] = lsm_flux[n];
330  }
331 
332  // Check if there is a user-specified roughness file to be read
333  std::string fname;
334  bool read_z0 = false;
335  if ( (flux_type == FluxCalcType::MOENG) ||
337  read_z0 = pp.query("most.roughness_file_name", fname);
338  }
339 
340  // Attributes for MFs and FABs
341  //--------------------------------------------------------
342  // Create a 2D ba, dm, & ghost cells
343  amrex::BoxArray ba = mf.boxArray();
344  amrex::BoxList bl2d = ba.boxList();
345  for (auto& b : bl2d) {
346  b.setRange(2, 0);
347  }
348  amrex::BoxArray ba2d(std::move(bl2d));
349  const amrex::DistributionMapping& dm = mf.DistributionMap();
350  const int ncomp = 1;
351  amrex::IntVect ng = mf.nGrowVect();
352  ng[2] = 0;
353 
354  // Z0 heights FAB
355  //--------------------------------------------------------
356  z_0[lev].define(ba2d, dm, ncomp, ng);
357  z_0[lev].setVal(z0_const);
358  if (read_z0) {
359  read_custom_roughness(lev, fname);
360  }
361 
362  // 2D MFs for U*, T*, T_surf
363  //--------------------------------------------------------
364  u_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
365  u_star[lev]->setVal(1.E34);
366 
367  w_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
368  w_star[lev]->setVal(1.E34);
369 
370  t_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
371  t_star[lev]->setVal(0.0); // default to neutral
372 
373  q_star[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
374  q_star[lev]->setVal(0.0); // default to dry
375 
376  olen[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
377  olen[lev]->setVal(1.E34);
378 
379  pblh[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
380  pblh[lev]->setVal(1.E34);
381 
382  t_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
383  t_surf[lev]->setVal(default_land_surf_temp);
384 
385  q_surf[lev] = std::make_unique<amrex::MultiFab>(ba2d, dm, ncomp, ng);
386  q_surf[lev]->setVal(default_land_surf_moist);
387 
388  // TODO: Do we want an enum struct for indexing?
389  if (m_sst_lev[lev][0] || m_tsk_lev[lev][0] || m_has_lsm_tsurf) {
390  // Valid SST, TSK or LSM data; t_surf set before computing fluxes (avoids
391  // extended lambda capture) Note that land temp will be set from m_tsk_lev
392  // while sea temp will be set from m_sst_lev
394  }
395 
396 
397  }
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:511
amrex::Vector< std::string > m_lsm_data_name
Definition: ERF_SurfaceLayer.H:612
bool m_has_lsm_tsurf
Definition: ERF_SurfaceLayer.H:594
void read_custom_roughness(const int &lev, const std::string &fname)
Definition: ERF_SurfaceLayer.cpp:676
@ 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 
)
678 {
679  // Read the file if we are on the coarsest level
680  if (lev==0) {
681  // Only the ioproc reads the file
682  Gpu::HostVector<Real> m_x,m_y,m_z0;
683  if (ParallelDescriptor::IOProcessor()) {
684  Print()<<"Reading MOST roughness file: "<< fname << std::endl;
685  std::ifstream file(fname);
686  Real value1,value2,value3;
687  while(file>>value1>>value2>>value3){
688  m_x.push_back(value1);
689  m_y.push_back(value2);
690  m_z0.push_back(value3);
691  }
692  file.close();
693  }
694 
695  // Broadcast the whole domain to every rank
696  int ioproc = ParallelDescriptor::IOProcessorNumber();
697  ParallelDescriptor::Barrier();
698  ParallelDescriptor::Bcast(m_x.data() , m_x.size() , ioproc);
699  ParallelDescriptor::Bcast(m_y.data() , m_x.size() , ioproc);
700  ParallelDescriptor::Bcast(m_z0.data(), m_z0.size(), ioproc);
701 
702  // Copy data to the GPU
703  int nnode = m_x.size();
704  Gpu::DeviceVector<Real> d_x(nnode),d_y(nnode),d_z0(nnode);
705  Gpu::copy(Gpu::hostToDevice, m_x.begin(), m_x.end(), d_x.begin());
706  Gpu::copy(Gpu::hostToDevice, m_y.begin(), m_y.end(), d_y.begin());
707  Gpu::copy(Gpu::hostToDevice, m_z0.begin(), m_z0.end(), d_z0.begin());
708  Real* xp = d_x.data();
709  Real* yp = d_y.data();
710  Real* z0p = d_z0.data();
711 
712  // Each rank populates it's z_0[lev] MultiFab
713  for (MFIter mfi(z_0[lev]); mfi.isValid(); ++mfi)
714  {
715  Box gtbx = mfi.growntilebox();
716 
717  // Populate z_phys data
718  Real tol = 1.0e-4;
719  auto dx = m_geom[lev].CellSizeArray();
720  auto ProbLoArr = m_geom[lev].ProbLoArray();
721  int ilo = m_geom[lev].Domain().smallEnd(0);
722  int jlo = m_geom[lev].Domain().smallEnd(1);
723  int klo = 0;
724  int ihi = m_geom[lev].Domain().bigEnd(0);
725  int jhi = m_geom[lev].Domain().bigEnd(1);
726 
727  Array4<Real> const& z0_arr = z_0[lev].array(mfi);
728  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
729  {
730  // Clip indices for ghost-cells
731  int ii = amrex::min(amrex::max(i,ilo),ihi);
732  int jj = amrex::min(amrex::max(j,jlo),jhi);
733 
734  // Location of nodes
735  Real x = ProbLoArr[0] + ii * dx[0];
736  Real y = ProbLoArr[1] + jj * dx[1];
737  int inode = ii + jj * (ihi-ilo+2); // stride is Nx+1
738  if (std::sqrt(std::pow(x-xp[inode],2)+std::pow(y-yp[inode],2)) < tol) {
739  z0_arr(i,j,klo) = z0p[inode];
740  } else {
741  // Unexpected list order, do brute force search
742  Real z0loc = 0.0;
743  bool found = false;
744  for (int n=0; n<nnode; ++n) {
745  Real delta=std::sqrt(std::pow(x-xp[n],2)+std::pow(y-yp[n],2));
746  if (delta < tol) {
747  found = true;
748  z0loc = z0p[n];
749  break;
750  }
751  }
752  AMREX_ASSERT_WITH_MESSAGE(found, "Location read from terrain file does not match the grid!");
753  amrex::ignore_unused(found);
754  z0_arr(i,j,klo) = z0loc;
755  }
756  });
757  } // mfi
758  } else {
759  // Create a BC mapper that uses FOEXTRAP at domain bndry
760  Vector<int> bc_lo(3,ERFBCType::foextrap);
761  Vector<int> bc_hi(3,ERFBCType::foextrap);
762  Vector<BCRec> bcr; bcr.push_back(BCRec(bc_lo.data(),bc_hi.data()));
763 
764  // Create ref ratio
765  IntVect ratio;
766  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
767  ratio[idim] = m_geom[lev].Domain().length(idim) / m_geom[0].Domain().length(idim);
768  }
769 
770  // Create interp object and interpolate from the coarsest grid
771  MFInterpolater* interp = &mf_cell_cons_interp;
772  interp->interp(z_0[0] , 0,
773  z_0[lev], 0,
774  1, z_0[lev].nGrowVect(),
775  m_geom[0], m_geom[lev],
776  m_geom[lev].Domain(),ratio,
777  bcr, 0);
778  }
779 }
@ 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,
const amrex::MultiFab &  cons_in,
const std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
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
17 {
18  // Update with SST/TSK data if we have a valid pointer
19  if (m_sst_lev[lev][0]) {
20  fill_tsurf_with_sst_and_tsk(lev, time);
21  }
22 
23  // Apply heating rate if needed
25  update_surf_temp(time);
26  }
27 
28  // Update qsurf with qsat over sea
30  fill_qsurf_with_qsat(lev, cons_in, z_phys_nd);
31  }
32 
33  // Update land surface temp if we have a valid pointer
34  if (m_has_lsm_tsurf) { get_lsm_tsurf(lev); }
35 
36  // Fill interior ghost cells
37  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
38 
39  // Compute plane averages for all vars (regardless of flux type)
41 
42  // Do we have a constant flux for moisture?
43  bool cons_qflux = ( (moist_type == MoistCalcType::MOISTURE_FLUX) ||
45 
46  // ***************************************************************
47  // Iterate the fluxes if moeng type
48  // First iterate over land -- the only model for surface roughness
49  // over land is RoughCalcType::CONSTANT
50  // ***************************************************************
53  bool is_land = true;
56  surface_flux most_flux(m_ma.get_zref(), surf_temp_flux, surf_moist_flux, cons_qflux);
57  compute_fluxes(lev, max_iters, most_flux, is_land);
58  } else {
59  amrex::Abort("Unknown value for rough_type_land");
60  }
63  surface_temp most_flux(m_ma.get_zref(), surf_temp_flux, surf_moist_flux, cons_qflux);
64  compute_fluxes(lev, max_iters, most_flux, is_land);
65  } else {
66  amrex::Abort("Unknown value for rough_type_land");
67  }
68  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
72  compute_fluxes(lev, max_iters, most_flux, is_land);
73  } else {
74  amrex::Abort("Unknown value for rough_type_land");
75  }
76  } else {
77  amrex::Abort("Unknown value for theta_type");
78  }
79  } // MOENG -- LAND
80 
81  // ***************************************************************
82  // Iterate the fluxes if moeng type
83  // Next iterate over sea -- the models for surface roughness
84  // over sea are CHARNOCK, DONELAN, MODIFIED_CHARNOCK or WAVE_COUPLED
85  // ***************************************************************
88  bool is_land = false;
93  cnk_a, cnk_visc, cons_qflux);
94  compute_fluxes(lev, max_iters, most_flux, is_land);
98  depth, cons_qflux);
99  compute_fluxes(lev, max_iters, most_flux, is_land);
100  } else if (rough_type_sea == RoughCalcType::DONELAN) {
101  surface_flux_donelan most_flux(m_ma.get_zref(),
103  cons_qflux);
104  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 
116  surface_temp_charnock most_flux(m_ma.get_zref(),
118  cnk_a, cnk_visc, cons_qflux);
119  compute_fluxes(lev, max_iters, most_flux, is_land);
123  depth, cons_qflux);
124  compute_fluxes(lev, max_iters, most_flux, is_land);
125  } else if (rough_type_sea == RoughCalcType::DONELAN) {
126  surface_temp_donelan most_flux(m_ma.get_zref(),
128  cons_qflux);
129  compute_fluxes(lev, max_iters, most_flux, is_land);
133  cons_qflux);
134  compute_fluxes(lev, max_iters, most_flux, is_land);
135  } else {
136  amrex::Abort("Unknown value for rough_type_sea");
137  }
138 
139  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
142  adiabatic_charnock most_flux(m_ma.get_zref(),
144  cnk_a, cnk_visc);
145  compute_fluxes(lev, max_iters, most_flux, is_land);
149  depth);
150  compute_fluxes(lev, max_iters, most_flux, is_land);
151  } else if (rough_type_sea == RoughCalcType::DONELAN) {
153  compute_fluxes(lev, max_iters, most_flux, is_land);
156  compute_fluxes(lev, max_iters, most_flux, is_land);
157  } else {
158  amrex::Abort("Unknown value for rough_type_sea");
159  }
160  } else {
161  amrex::Abort("Unknown value for theta_type");
162  }
163 
164  } // MOENG -- SEA
165 
167  if (custom_rhosurf > 0) {
168  specified_rho_surf = true;
169  u_star[lev]->setVal(std::sqrt(custom_rhosurf) * custom_ustar);
170  t_star[lev]->setVal(custom_rhosurf * custom_tstar);
171  q_star[lev]->setVal(custom_rhosurf * custom_qstar);
172  } else {
173  u_star[lev]->setVal(custom_ustar);
174  t_star[lev]->setVal(custom_tstar);
175  q_star[lev]->setVal(custom_qstar);
176  }
177  }
178 }
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:679
void update_surf_temp(const amrex::Real &time)
Definition: ERF_SurfaceLayer.H:461
void get_lsm_tsurf(const int &lev)
Definition: ERF_SurfaceLayer.cpp:617
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:577
void compute_fluxes(const int &lev, const int &max_iters, const FluxIter &most_flux, bool is_land)
Definition: ERF_SurfaceLayer.cpp:189
void fill_tsurf_with_sst_and_tsk(const int &lev, const amrex::Real &time)
Definition: ERF_SurfaceLayer.cpp:505
Definition: ERF_MOSTStress.H:187
Definition: ERF_MOSTStress.H:357
Definition: ERF_MOSTStress.H:276
Definition: ERF_MOSTStress.H:427
Definition: ERF_MOSTStress.H:137
Definition: ERF_MOSTStress.H:608
Definition: ERF_MOSTStress.H:838
Definition: ERF_MOSTStress.H:727
Definition: ERF_MOSTStress.H:946
Definition: ERF_MOSTStress.H:503
Definition: ERF_MOSTStress.H:1197
Definition: ERF_MOSTStress.H:1521
Definition: ERF_MOSTStress.H:1363
Definition: ERF_MOSTStress.H:1639
Definition: ERF_MOSTStress.H:1063

◆ 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
478  {
479  m_ma.update_field_ptrs(lev, vars_old, Theta_prim, Qv_prim, Qr_prim);
480  }
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 MoistureComponentIndices moisture_indices 
)
653 {
655  MYNNPBLH estimator;
656  compute_pblh(lev, vars, z_phys_cc, estimator, moisture_indices);
658  amrex::Error("YSU/MRF PBLH calc not implemented yet");
659  }
660 }
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_surf_temp()

void SurfaceLayer::update_surf_temp ( const amrex::Real time)
inline
462  {
463  if (surf_heating_rate != 0) {
464  int nlevs = static_cast<int>(m_geom.size());
465  for (int lev = 0; lev < nlevs; lev++) {
466  t_surf[lev]->setVal(surf_temp + surf_heating_rate * time);
467  amrex::Print() << "Surface temp at t=" << time << ": "
468  << surf_temp + surf_heating_rate * time << std::endl;
469  }
470  }
471  }

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}

◆ 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

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

◆ 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_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: