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

#include <ERF_SurfaceLayer.H>

Collaboration diagram for SurfaceLayer:

Public Types

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

Public Member Functions

 SurfaceLayer (const amrex::Vector< amrex::Geometry > &geom, bool &use_rot_surface_flux, std::string a_pp_prefix, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &z_phys_nd, const MeshType &a_mesh_type, const TerrainType &a_terrain_type, const TurbChoice &a_turb_choice, amrex::Real start_low_time, amrex::Real final_low_time, amrex::Real low_time_interval=zero, const amrex::Vector< const eb_ * > &eb_vec={})
 
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 &elapsed_time, const amrex::Real &elapsed_time_since_start_low, amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, const std::unique_ptr< amrex::MultiFab > &walldist, int max_iters=100)
 
template<typename FluxIter >
void compute_fluxes (const int &lev, const int &max_iters, amrex::MultiFab &cons_in, const FluxIter &most_flux, bool is_land)
 
void init_tke_from_ustar (const int &lev, amrex::MultiFab &cons, const std::unique_ptr< amrex::MultiFab > &z_phys_nd, const amrex::Real tkefac=one, const amrex::Real zscale=amrex::Real(700.0))
 
void impose_SurfaceLayer_bcs (const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Tau_lev, amrex::MultiFab *xheat_flux, amrex::MultiFab *yheat_flux, amrex::MultiFab *zheat_flux, amrex::MultiFab *xqv_flux, amrex::MultiFab *yqv_flux, amrex::MultiFab *zqv_flux, const amrex::MultiFab *z_phys)
 
void impose_SurfaceLayer_bcs_EB (const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< 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)
 
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)
 
template<typename FluxCalc >
void compute_SurfaceLayer_bcs_EB (const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< 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 FluxCalc &flux_comp)
 
void compute_sfc_params_from_lsm_fluxes (const int &lev, amrex::MultiFab &cons_in)
 
void fill_tsurf_with_sst_and_tsk (const int &lev, const amrex::Real &time)
 
void fill_qsurf_with_qsat (const int &lev, const amrex::MultiFab &cons_in, const std::unique_ptr< amrex::MultiFab > &z_phys_nd)
 
void get_lsm_tsurf (const int &lev)
 
void update_pblh (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const MoistureComponentIndices &moisture_indices)
 
template<typename PBLHeightEstimator >
void compute_pblh (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indice)
 
void read_custom_roughness (const int &lev, const std::string &fname)
 
void update_surf_temp (const amrex::Real &time)
 
void update_mac_ptrs (const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars_old, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Theta_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qv_prim, amrex::Vector< std::unique_ptr< amrex::MultiFab >> &Qr_prim)
 
amrex::MultiFab * get_u_star (const int &lev)
 
amrex::MultiFab * get_w_star (const int &lev)
 
amrex::MultiFab * get_t_star (const int &lev)
 
amrex::MultiFab * get_q_star (const int &lev)
 
amrex::MultiFab * get_olen (const int &lev)
 
amrex::MultiFab * get_pblh (const int &lev)
 
const amrex::MultiFab * get_mac_avg (const int &lev, int comp)
 
amrex::MultiFab * get_t_surf (const int &lev)
 
void set_t_surf (const int &lev, const amrex::Real tsurf)
 
amrex::MultiFab * get_q_surf (const int &lev)
 
void set_q_surf (const int &lev, const amrex::Real qsurf)
 
amrex::Real get_zref (const int &lev)
 
amrex::MultiFab * get_z0 (const int &lev)
 
bool have_variable_sea_roughness ()
 
amrex::iMultiFab * get_lmask (const int &lev)
 
int lmask_min_reduce (amrex::iMultiFab &lmask, const int &nghost)
 
void update_sst_ptr (const int lev, const int itime, amrex::MultiFab *sst_ptr)
 
void update_tsk_ptr (const int lev, const int itime, amrex::MultiFab *tsk_ptr)
 
template<typename FluxIter >
void compute_fluxes (const int &lev, const int &max_iters, MultiFab &cons_in, const FluxIter &most_flux, bool is_land)
 
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 FluxCalc >
void compute_SurfaceLayer_bcs_EB (const int &lev, Vector< const MultiFab * > mfs, Vector< Vector< std::unique_ptr< MultiFab >>> &Tau_EB, [[maybe_unused]] MultiFab *xheat_flux, [[maybe_unused]] MultiFab *yheat_flux, MultiFab *Hfx3_EB, [[maybe_unused]] MultiFab *xqv_flux, [[maybe_unused]] MultiFab *yqv_flux, [[maybe_unused]] MultiFab *zqv_flux, const FluxCalc &flux_comp)
 
template<typename PBLHeightEstimator >
void compute_pblh (const int &lev, Vector< Vector< MultiFab >> &vars, MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indices)
 

Public Attributes

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

Private Attributes

amrex::Vector< amrex::Geometry > m_geom
 
bool m_rotate = false
 
amrex::Real m_start_low_time
 
amrex::Real m_final_low_time
 
amrex::Real m_low_time_interval
 
bool m_include_wstar = false
 
amrex::Real z0_const {amrex::Real(0.1)}
 
amrex::Real default_land_surf_temp {amrex::Real(300.)}
 
amrex::Real surf_temp
 
amrex::Real surf_heating_rate {0}
 
amrex::Real surf_temp_flux {0}
 
amrex::Real default_land_surf_moist {zero}
 
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 {amrex::Real(0.0185)}
 
bool cnk_visc {false}
 
amrex::Real depth {amrex::Real(30.0)}
 
amrex::Vector< amrex::MultiFab > z_0
 
bool m_var_z0 {false}
 
amrex::Real rico_theta_z0 {amrex::Real(298.0)}
 
amrex::Real rico_qsat_z0 {amrex::Real(0.001)}
 
bool use_moisture
 
bool m_has_lsm_fluxes = false
 
bool m_has_lsm_tsurf = false
 
bool m_has_ocean_lsm_tsurf = false
 
int m_lsm_tsurf_indx = -1
 
amrex::Real m_Cd = zero
 
amrex::Real m_Ch = zero
 
amrex::Real m_Cq = zero
 
bool m_ignore_sst = false
 
amrex::Vector< const eb_ * > m_eb_vec
 
TerrainType m_terrain_type
 
MOSTAverage m_ma
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
 
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
 
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
 
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
 
amrex::Vector< std::string > m_lsm_data_name
 
amrex::Vector< std::string > m_lsm_flux_name
 
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
 
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
 
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
 
bool m_update_k_rans = false
 
amrex::Real inv_Cmu2 = zero
 
amrex::Real theta_ref = zero
 

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–amrex::Real(489.) https://doi.org/amrex::Real(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/amrex::Real(10.1142)/amrex::Real(2975.)

Member Enumeration Documentation

◆ FluxCalcType

Enumerator
MOENG 

Moeng functional form.

DONELAN 

Donelan functional form.

CUSTOM 

Custom constant flux functional form.

BULK_COEFF 

Bulk transfer coefficient functional form.

ROTATE 

Terrain rotation flux functional form.

RICO 
678  {
679  MOENG = 0, ///< Moeng functional form
680  DONELAN, ///< Donelan functional form
681  CUSTOM, ///< Custom constant flux functional form
682  BULK_COEFF, ///< Bulk transfer coefficient functional form
683  ROTATE, ///< Terrain rotation flux functional form
684  RICO
685  };

◆ MoistCalcType

Enumerator
ADIABATIC 
MOISTURE_FLUX 

Qv-flux specified.

SURFACE_MOISTURE 

Surface Qv specified.

693  {
694  ADIABATIC = 0,
695  MOISTURE_FLUX, ///< Qv-flux specified
696  SURFACE_MOISTURE ///< Surface Qv specified
697  };

◆ PBLHeightCalcType

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

◆ RoughCalcType

Enumerator
CONSTANT 

Constant z0.

CHARNOCK 
MODIFIED_CHARNOCK 
DONELAN 
WAVE_COUPLED 
699  {
700  CONSTANT = 0, ///< Constant z0
701  CHARNOCK,
702  MODIFIED_CHARNOCK,
703  DONELAN,
704  WAVE_COUPLED
705  };

◆ ThetaCalcType

Enumerator
ADIABATIC 
HEAT_FLUX 

Heat-flux specified.

SURFACE_TEMPERATURE 

Surface temperature specified.

687  {
688  ADIABATIC = 0,
689  HEAT_FLUX, ///< Heat-flux specified
690  SURFACE_TEMPERATURE ///< Surface temperature specified
691  };

Constructor & Destructor Documentation

◆ SurfaceLayer()

SurfaceLayer::SurfaceLayer ( const amrex::Vector< amrex::Geometry > &  geom,
bool &  use_rot_surface_flux,
std::string  a_pp_prefix,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  Qv_prim,
amrex::Vector< std::unique_ptr< amrex::MultiFab >> &  z_phys_nd,
const MeshType &  a_mesh_type,
const TerrainType &  a_terrain_type,
const TurbChoice a_turb_choice,
amrex::Real  start_low_time,
amrex::Real  final_low_time,
amrex::Real  low_time_interval = zero,
const amrex::Vector< const eb_ * > &  eb_vec = {} 
)
inlineexplicit
48  {})
49  : m_geom(geom),
50  m_rotate(use_rot_surface_flux),
51  m_start_low_time(start_low_time),
52  m_final_low_time(final_low_time),
53  m_low_time_interval(low_time_interval),
54  m_eb_vec(eb_vec),
55  m_terrain_type(a_terrain_type),
56  m_ma(geom, (z_phys_nd[0] != nullptr), a_pp_prefix, a_mesh_type, a_terrain_type, eb_vec)
57  {
58  // We have a moisture model if Qv_prim is a valid pointer
59  use_moisture = (Qv_prim[0].get());
60 
61  // Get roughness
62  amrex::ParmParse pp("erf");
63  pp.query("most.z0", z0_const);
64 
65  // Specify how to compute the flux
66  if (use_rot_surface_flux) {
68  } else {
69  std::string flux_string_in;
70  std::string flux_string{"moeng"};
71  auto read_flux = pp.query("surface_layer.flux_type", flux_string_in);
72  if (read_flux) {
73  flux_string = amrex::toLower(flux_string_in);
74  }
75  if (flux_string == "donelan") {
77  } else if (flux_string == "moeng") {
79  } else if (flux_string == "rico") {
81  } else if (flux_string == "bulk_coeff") {
83  } else if (flux_string == "custom") {
85  } else {
86  amrex::Abort("Undefined MOST flux type!");
87  }
88  }
89 
90  // Include w* to handle free convection (Beljaars 1995, QJRMS)
91  pp.query("most.include_wstar", m_include_wstar);
92 
93  std::string pblh_string_in;
94  std::string pblh_string{"none"};
95  auto read_pblh = pp.query("most.pblh_calc", pblh_string_in);
96  if (read_pblh) {
97  pblh_string = amrex::toLower(pblh_string_in);
98  }
99  if (pblh_string == "none") {
101  } else if (pblh_string == "mynn25") {
103  } else if (pblh_string == "mynnedmf") {
105  } else if (pblh_string == "ysu") {
107  } else if (pblh_string == "mrf") {
109  } else {
110  amrex::Abort("Undefined PBLH calc type!");
111  }
112 
113  // Get surface temperature
114  auto erf_st = pp.query("most.surf_temp", surf_temp);
115  if (erf_st) { default_land_surf_temp = surf_temp; }
116 
117  // Get surface moisture
118  bool erf_sq = false;
119  if (use_moisture) { erf_sq = pp.query("most.surf_moist", surf_moist); }
120  if (erf_sq) { default_land_surf_moist = surf_moist; }
121 
122  // Custom type user must specify the fluxes
127  pp.get("most.ustar", custom_ustar);
128  pp.get("most.tstar", custom_tstar);
129  pp.get("most.qstar", custom_qstar);
130  pp.query("most.rhosurf", custom_rhosurf);
131  if (custom_qstar != 0) {
133  "Specified custom MOST qv flux without moisture model!");
134  }
135  amrex::Print() << "Using specified ustar, tstar, qstar for MOST = "
136  << custom_ustar << " " << custom_tstar << " "
137  << custom_qstar << std::endl;
138 
139  // Bulk transfer coefficient (must specify coeffs and surface values)
140  } else if (flux_type == FluxCalcType::BULK_COEFF) {
141  pp.get("most.Cd", m_Cd);
142  pp.get("most.Ch", m_Ch);
143  pp.get("most.Cq", m_Cq);
144  pp.get("most.surf_temp", default_land_surf_temp);
145  pp.get("most.surf_moist", default_land_surf_moist);
146  amrex::Print() << "Using specified Cd, Ch, Cq for MOST = "
147  << m_Cd << " " << m_Ch << " "
148  << m_Cq << std::endl;
149 
150  // Specify surface temperature/moisture or surface flux
151  } else {
152  if (erf_st) {
154  pp.query("most.surf_heating_rate", surf_heating_rate); // [K/h]
155 
156  // Modify rate to be in units of K / s rather than K / hr
157  surf_heating_rate /= amrex::Real(3600.0); // [K/s]
158 
159  if (pp.query("most.surf_temp_flux", surf_temp_flux)) {
160  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
161  }
162  } else {
163  pp.query("most.surf_temp_flux", surf_temp_flux);
164 
165  if (pp.query("most.surf_heating_rate", surf_heating_rate)) {
166  amrex::Abort("Can only specify one of surf_temp_flux or surf_heating_rate");
167  }
168  if (std::abs(surf_temp_flux) >
171  } else {
173  }
174  }
175 
176  if (erf_sq) {
178  } else {
179  pp.query("most.surf_moist_flux", surf_moist_flux);
180  if (std::abs(surf_moist_flux) >
183  } else {
185  }
186  }
187  }
188 
190  {
191  pp.query("most.rico.theta_z0", rico_theta_z0);
192  pp.query("most.rico.qsat_z0", rico_qsat_z0);
193  }
194 
195  // Make sure the inputs file doesn't try to use most.roughness_type
196  std::string bogus_input;
197  if (pp.query("most.roughness_type", bogus_input) > 0) {
198  amrex::Abort("most.roughness_type is deprecated; use "
199  "most.roughness_type_land and/or most.roughness_type_sea");
200  }
201 
202  // Specify how to compute the surface flux over land (if there is any)
203  std::string rough_land_string_in;
204  std::string rough_land_string{"constant"};
205  auto read_rough_land =
206  pp.query("most.roughness_type_land", rough_land_string_in);
207  if (read_rough_land) {
208  rough_land_string = amrex::toLower(rough_land_string_in);
209  }
210  if (rough_land_string == "constant") {
212  } else {
213  amrex::Abort("Undefined MOST roughness type for land!");
214  }
215 
216  // Specify how to compute the surface flux over sea (if there is any)
217  std::string rough_sea_string_in;
218  std::string rough_sea_string{"charnock"};
219  auto read_rough_sea = pp.query("most.roughness_type_sea", rough_sea_string_in);
220  if (read_rough_sea) {
221  rough_sea_string = amrex::toLower(rough_sea_string_in);
222  }
223  if (rough_sea_string == "charnock") {
225  pp.query("most.charnock_constant", cnk_a);
226  pp.query("most.charnock_viscosity", cnk_visc);
227  if (cnk_a > 0) {
228  amrex::Print() << "If there is water, Charnock relation with C_a="
229  << cnk_a << (cnk_visc ? " and viscosity" : "")
230  << " will be used" << std::endl;
231  } else {
232  amrex::Print() << "If there is water, Charnock relation with variable "
233  "Charnock parameter (COARE3.0)"
234  << (cnk_visc ? " and viscosity" : "") << " will be used"
235  << std::endl;
236  }
237  } else if (rough_sea_string == "coare3.0") {
239  amrex::Print() << "If there is water, Charnock relation with variable "
240  "Charnock parameter (COARE3.0)"
241  << (cnk_visc ? " and viscosity" : "") << " will be used"
242  << std::endl;
243  cnk_a = -1;
244  } else if (rough_sea_string == "donelan") {
246  } else if (rough_sea_string == "modified_charnock") {
248  pp.query("most.modified_charnock_depth", depth);
249  } else if (rough_sea_string == "wave_coupled") {
251  } else if (rough_sea_string == "constant") {
253  } else {
254  amrex::Abort("Undefined MOST roughness type for sea!");
255  }
256 
257  // use skin temperature instead of sea-surface temperature
258  // (wrfinput data may have lower resolution SST data)
259  pp.query("most.ignore_sst", m_ignore_sst);
260 
261  // If we're using the RANS k model, then we need to update the dirichlet
262  // BC based on the instantaneous u* and θ*; the turbulence modeling
263  // choices can vary per level but for now, assume that if specified then
264  // all levels are using the same RANS model.
265  m_update_k_rans = (a_turb_choice.rans_type == RANSType::kEqn &&
266  a_turb_choice.dirichlet_k == true);
267  if (m_update_k_rans) {
268  inv_Cmu2 = one / (a_turb_choice.Cmu0 * a_turb_choice.Cmu0);
269  theta_ref = a_turb_choice.theta_ref;
270  }
271 
272  } // constructor
constexpr amrex::Real one
Definition: ERF_Constants.H:7
ParmParse pp("prob")
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_ASSERT_WITH_MESSAGE(wbar_cutoff_min > wbar_cutoff_max, "ERROR: wbar_cutoff_min < wbar_cutoff_max")
ThetaCalcType theta_type
Definition: ERF_SurfaceLayer.H:710
bool m_include_wstar
Definition: ERF_SurfaceLayer.H:724
bool m_rotate
Definition: ERF_SurfaceLayer.H:719
PBLHeightCalcType pblh_type
Definition: ERF_SurfaceLayer.H:714
bool use_moisture
Definition: ERF_SurfaceLayer.H:747
amrex::Real m_Cq
Definition: ERF_SurfaceLayer.H:755
amrex::Vector< const eb_ * > m_eb_vec
Definition: ERF_SurfaceLayer.H:759
RoughCalcType rough_type_land
Definition: ERF_SurfaceLayer.H:712
amrex::Real z0_const
Definition: ERF_SurfaceLayer.H:725
amrex::Real cnk_a
Definition: ERF_SurfaceLayer.H:738
amrex::Real m_Ch
Definition: ERF_SurfaceLayer.H:754
amrex::Real m_start_low_time
Definition: ERF_SurfaceLayer.H:720
amrex::Real surf_temp
Definition: ERF_SurfaceLayer.H:727
amrex::Real rico_qsat_z0
Definition: ERF_SurfaceLayer.H:745
bool m_update_k_rans
Definition: ERF_SurfaceLayer.H:782
amrex::Real surf_moist_flux
Definition: ERF_SurfaceLayer.H:732
RoughCalcType rough_type_sea
Definition: ERF_SurfaceLayer.H:713
amrex::Real surf_moist
Definition: ERF_SurfaceLayer.H:731
bool m_ignore_sst
Definition: ERF_SurfaceLayer.H:757
amrex::Real custom_qstar
Definition: ERF_SurfaceLayer.H:735
amrex::Real custom_rhosurf
Definition: ERF_SurfaceLayer.H:736
@ MOENG
Moeng functional form.
@ BULK_COEFF
Bulk transfer coefficient functional form.
@ CUSTOM
Custom constant flux functional form.
@ ROTATE
Terrain rotation flux functional form.
@ DONELAN
Donelan functional form.
@ SURFACE_MOISTURE
Surface Qv specified.
@ MOISTURE_FLUX
Qv-flux specified.
amrex::Real depth
Definition: ERF_SurfaceLayer.H:740
amrex::Real default_land_surf_moist
Definition: ERF_SurfaceLayer.H:730
amrex::Real rico_theta_z0
Definition: ERF_SurfaceLayer.H:744
amrex::Real surf_temp_flux
Definition: ERF_SurfaceLayer.H:729
amrex::Vector< amrex::Geometry > m_geom
Definition: ERF_SurfaceLayer.H:718
amrex::Real theta_ref
Definition: ERF_SurfaceLayer.H:784
amrex::Real custom_tstar
Definition: ERF_SurfaceLayer.H:734
amrex::Real m_final_low_time
Definition: ERF_SurfaceLayer.H:721
bool cnk_visc
Definition: ERF_SurfaceLayer.H:739
amrex::Real surf_heating_rate
Definition: ERF_SurfaceLayer.H:728
FluxCalcType flux_type
Definition: ERF_SurfaceLayer.H:709
MoistCalcType moist_type
Definition: ERF_SurfaceLayer.H:711
amrex::Real inv_Cmu2
Definition: ERF_SurfaceLayer.H:783
amrex::Real custom_ustar
Definition: ERF_SurfaceLayer.H:733
amrex::Real m_Cd
Definition: ERF_SurfaceLayer.H:753
amrex::Real default_land_surf_temp
Definition: ERF_SurfaceLayer.H:726
@ SURFACE_TEMPERATURE
Surface temperature specified.
@ HEAT_FLUX
Heat-flux specified.
TerrainType m_terrain_type
Definition: ERF_SurfaceLayer.H:760
amrex::Real m_low_time_interval
Definition: ERF_SurfaceLayer.H:722
MOSTAverage m_ma
Definition: ERF_SurfaceLayer.H:761
real(c_double), parameter epsilon
Definition: ERF_module_model_constants.F90:12
RANSType rans_type
Definition: ERF_TurbStruct.H:454
bool dirichlet_k
Definition: ERF_TurbStruct.H:456
amrex::Real Cmu0
Definition: ERF_TurbStruct.H:432
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:443

Member Function Documentation

◆ compute_fluxes() [1/2]

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

◆ compute_fluxes() [2/2]

template<typename FluxIter >
void SurfaceLayer::compute_fluxes ( const int &  lev,
const int &  max_iters,
MultiFab &  cons_in,
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
277 {
278  // Pointers to the computed averages
279  const auto *const tm_ptr = m_ma.get_average(lev,2); // potential temperature
280  const auto *const qvm_ptr = m_ma.get_average(lev,3); // water vapor mixing ratio
281  const auto *const tvm_ptr = m_ma.get_average(lev,4); // virtual potential temperature
282  const auto *const umm_ptr = m_ma.get_average(lev,5); // horizontal velocity magnitude
283  const auto *const zref_ptr = m_ma.get_zref(lev); // reference height
284  const bool l_use_eb = (m_terrain_type == TerrainType::EB);
285 
286  const int klo = m_geom[lev].Domain().smallEnd(2);
287  IntVect ng = u_star[lev]->nGrowVect(); ng[2] = 0;
288 
289  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi)
290  {
291  Box gtbx = mfi.tilebox(IntVect(0),ng);
292 
293  if (!l_use_eb && gtbx.smallEnd(2) != klo) { continue; }
294 
295  if (!l_use_eb) { gtbx.makeSlab(2,klo); }
296 
297  auto u_star_arr = u_star[lev]->array(mfi);
298  auto t_star_arr = t_star[lev]->array(mfi);
299  auto q_star_arr = q_star[lev]->array(mfi);
300  auto t_surf_arr = t_surf[lev]->array(mfi);
301  auto q_surf_arr = q_surf[lev]->array(mfi);
302  auto olen_arr = olen[lev]->array(mfi);
303 
304  const auto tm_arr = tm_ptr->array(mfi);
305  const auto tvm_arr = tvm_ptr->array(mfi);
306  const auto qvm_arr = qvm_ptr->array(mfi);
307  const auto umm_arr = umm_ptr->array(mfi);
308  const auto zref_arr = zref_ptr->array(mfi);
309  const auto z0_arr = z_0[lev].array(mfi);
310 
311  // PBL height if we need to calculate wstar for the Beljaars correction
312  // TODO: can/should we apply this in LES mode?
313  const auto w_star_arr = (m_include_wstar) ? w_star[lev].get()->array(mfi) : Array4<Real> {};
314  const auto pblh_arr = (m_include_wstar) ? pblh[lev].get()->array(mfi) : Array4<Real> {};
315 
316  // Wave properties if they exist
317  const auto Hwave_arr = (m_Hwave_lev[lev]) ? m_Hwave_lev[lev]->array(mfi) : Array4<Real> {};
318  const auto Lwave_arr = (m_Lwave_lev[lev]) ? m_Lwave_lev[lev]->array(mfi) : Array4<Real> {};
319  const auto eta_arr = (m_eddyDiffs_lev[lev]) ? m_eddyDiffs_lev[lev]->array(mfi) : Array4<Real> {};
320 
321  // Land mask array if it exists
322  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
323  Array4<int> {};
324 
325  // Get EB flags if needed
326  const auto flag_arr = (l_use_eb) ? m_eb_vec[lev]->get_const_factory()->getMultiEBCellFlagFab()[mfi].const_array() : Array4<const EBCellFlag>{};
327 
328  if (!l_use_eb) {
329  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
330  {
331  if (( is_land && lmask_arr(i,j,0) == 1) ||
332  (!is_land && lmask_arr(i,j,0) == 0))
333  {
334  // NOTE: All 2D MFs so k index is always 0 from ba2d definition
335  most_flux.iterate_flux(i, j, 0, max_iters,
336  zref_arr, // set in most average
337  z0_arr, // updated if(!is_land)
338  umm_arr, tm_arr, tvm_arr, qvm_arr,
339  u_star_arr, // updated
340  w_star_arr, // updated if(m_include_wstar)
341  t_star_arr, q_star_arr, // updated
342  t_surf_arr, q_surf_arr, olen_arr, // updated
343  pblh_arr, // updated if(m_include_wstar)
344  Hwave_arr, Lwave_arr, eta_arr);
345  }
346  });
347  // EB
348  } else {
352  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
353  {
354  if (( is_land && lmask_arr(i,j,0) == 1) ||
355  (!is_land && lmask_arr(i,j,0) == 0))
356  {
357  if (flag_arr(i,j,k).isSingleValued()) {
358  most_flux.iterate_flux(i, j, k, max_iters,
359  zref_arr, // set in most average
360  z0_arr, // updated if(!is_land)
361  umm_arr, tm_arr, tvm_arr, qvm_arr,
362  u_star_arr, // updated
363  w_star_arr, // updated if(m_include_wstar)
364  t_star_arr, q_star_arr, // updated
365  t_surf_arr, q_surf_arr, olen_arr, // updated
366  pblh_arr, // updated if(m_include_wstar)
367  Hwave_arr, Lwave_arr, eta_arr);
368  }
369  }
370  });
371  } else {
372  amrex::Abort("FluxIter type not supported for EB");
373  }
374  }
375  }
376 }
amrex::Real value
Definition: ERF_HurricaneDiagnostics.H:20
pp get("wavelength", wavelength)
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
amrex::MultiFab * get_zref(const int &lev) const
Definition: ERF_MOSTAverage.H:116
const amrex::MultiFab * get_average(const int &lev, const int &comp) const
Definition: ERF_MOSTAverage.H:113
amrex::Vector< amrex::Vector< amrex::iMultiFab * > > m_lmask_lev
Definition: ERF_SurfaceLayer.H:773
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_surf
Definition: ERF_SurfaceLayer.H:768
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_star
Definition: ERF_SurfaceLayer.H:765
amrex::Vector< amrex::MultiFab * > m_Lwave_lev
Definition: ERF_SurfaceLayer.H:779
amrex::Vector< amrex::MultiFab > z_0
Definition: ERF_SurfaceLayer.H:741
amrex::Vector< std::unique_ptr< amrex::MultiFab > > w_star
Definition: ERF_SurfaceLayer.H:763
amrex::Vector< std::unique_ptr< amrex::MultiFab > > u_star
Definition: ERF_SurfaceLayer.H:762
amrex::Vector< std::unique_ptr< amrex::MultiFab > > q_surf
Definition: ERF_SurfaceLayer.H:769
amrex::Vector< amrex::MultiFab * > m_Hwave_lev
Definition: ERF_SurfaceLayer.H:778
amrex::Vector< std::unique_ptr< amrex::MultiFab > > t_star
Definition: ERF_SurfaceLayer.H:764
amrex::Vector< amrex::MultiFab * > m_eddyDiffs_lev
Definition: ERF_SurfaceLayer.H:780
amrex::Vector< std::unique_ptr< amrex::MultiFab > > olen
Definition: ERF_SurfaceLayer.H:766
amrex::Vector< std::unique_ptr< amrex::MultiFab > > pblh
Definition: ERF_SurfaceLayer.H:767
@ ng
Definition: ERF_Morrison.H:48
Here is the call graph for this function:

◆ 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 
)
1164 {
1165  est.compute_pblh(m_geom[lev],z_phys_cc, pblh[lev].get(),
1166  vars[lev][Vars::cons],m_lmask_lev[lev][0],
1167  moisture_indices);
1168 }
@ cons
Definition: ERF_IndexDefines.H:174
Here is the call graph for this function:

◆ compute_sfc_params_from_lsm_fluxes()

void SurfaceLayer::compute_sfc_params_from_lsm_fluxes ( const int &  lev,
amrex::MultiFab &  cons_in 
)
909 {
911  bool has_moisture = use_moisture;
912  const int klo = m_geom[lev].Domain().smallEnd(2);
913  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi) {
914 
915  Box vbx = mfi.validbox();
916  if (vbx.smallEnd(2) != klo) { continue; }
917  vbx.makeSlab(2,0);
918 
919  // Get CC state
920  const Array4<const Real> cons_arr = cons_in.const_array(mfi);
921 
922  // Get SL params
923  const auto u_star_arr = u_star[lev]->array(mfi);
924  const auto t_star_arr = t_star[lev]->array(mfi);
925  const auto q_star_arr = q_star[lev]->array(mfi);
926  const auto olen_arr = olen[lev]->array(mfi);
927 
928  // Get LSM fluxes
929  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
930  Array4<int> {};
931  auto lsm_t_flux_arr = Array4<Real> {};
932  auto lsm_q_flux_arr = Array4<Real> {};
933  auto lsm_tau13_arr = Array4<Real> {};
934  auto lsm_tau23_arr = Array4<Real> {};
935  for (int n(0); n<m_lsm_flux_lev[lev].size(); ++n) {
936  if (toLower(m_lsm_flux_name[n]) == "t_flux") { lsm_t_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
937  if (toLower(m_lsm_flux_name[n]) == "q_flux") { lsm_q_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
938  if (toLower(m_lsm_flux_name[n]) == "tau13") { lsm_tau13_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
939  if (toLower(m_lsm_flux_name[n]) == "tau23") { lsm_tau23_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
940  }
941 
942  ParallelFor(vbx, [=] AMREX_GPU_DEVICE(int i, int j, int /*k*/) noexcept
943  {
944  int is_land = (lmask_arr) ? lmask_arr(i,j,0) : 1;
945  // Skip cells the LSM did not flux (sea-ice / open water -> the LSM
946  // wrote the lsm_flux_undefined sentinel). They keep the u*/T*/q*/L
947  // from the MOST iteration (computed from the surface temperature).
948  // Using the sentinel here would give u_star ~ sqrt(1e30) and blow up
949  // the PBL on the next step.
950  if (is_land && lsm_t_flux_arr && lsm_t_flux_arr(i,j,0) < lsm_flux_undefined) {
951  Real rho = cons_arr(i,j,klo,Rho_comp);
952  Real Thd = cons_arr(i,j,klo,RhoTheta_comp) / rho;
953  Real qv = (has_moisture) ? cons_arr(i,j,klo,RhoQ1_comp) / rho : zero;
954  Real Thv = Thd * (one + (R_v/R_d - one)*qv);
955  Real tau = std::sqrt( lsm_tau13_arr(i,j,0)*lsm_tau13_arr(i,j,0)
956  + lsm_tau23_arr(i,j,0)*lsm_tau23_arr(i,j,0) );
957  u_star_arr(i,j,0) = amrex::max(std::sqrt(tau),eps);
958  if (lsm_t_flux_arr(i,j,0)>=zero) {
959  t_star_arr(i,j,0) = amrex::min(-lsm_t_flux_arr(i,j,0) / u_star_arr(i,j,0),-eps);
960  } else {
961  t_star_arr(i,j,0) = amrex::max(-lsm_t_flux_arr(i,j,0) / u_star_arr(i,j,0),eps);
962  }
963  if (lsm_q_flux_arr(i,j,0)>=zero) {
964  q_star_arr(i,j,0) = amrex::min(-lsm_q_flux_arr(i,j,0) / u_star_arr(i,j,0),-eps);
965  } else {
966  q_star_arr(i,j,0) = amrex::max(-lsm_q_flux_arr(i,j,0) / u_star_arr(i,j,0),eps);
967  }
968  olen_arr(i,j,0) = ( u_star_arr(i,j,0) * u_star_arr(i,j,0) * Thv ) /
969  ( KAPPA * CONST_GRAV * t_star_arr(i,j,0) );
970  }
971  });
972  } // mfi
973 }
constexpr amrex::Real R_v
Definition: ERF_Constants.H:30
constexpr amrex::Real KAPPA
Definition: ERF_Constants.H:39
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real lsm_flux_undefined
Definition: ERF_Constants.H:22
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:40
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
rho
Definition: ERF_InitCustomPert_Bubble.H:106
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_flux_lev
Definition: ERF_SurfaceLayer.H:775
amrex::Vector< std::string > m_lsm_flux_name
Definition: ERF_SurfaceLayer.H:777
@ qv
Definition: ERF_Kessler.H:29
Here is the call graph for this function:

◆ 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
489 {
490  bool rotate = m_rotate;
491  const int klo = m_geom[lev].Domain().smallEnd(2);
492  const auto& dxInv = m_geom[lev].InvCellSizeArray();
493  for (MFIter mfi(*mfs[0]); mfi.isValid(); ++mfi)
494  {
495  // Get field arrays
496  const auto cons_arr = mfs[Vars::cons]->array(mfi);
497  const auto velx_arr = mfs[Vars::xvel]->array(mfi);
498  const auto vely_arr = mfs[Vars::yvel]->array(mfi);
499  const auto velz_arr = mfs[Vars::zvel]->array(mfi);
500 
501  // Diffusive stress vars
502  auto t13_arr = Tau_lev[TauType::tau13]->array(mfi);
503  auto t31_arr = (Tau_lev[TauType::tau31]) ? Tau_lev[TauType::tau31]->array(mfi) : Array4<Real>{};
504 
505  auto t23_arr = Tau_lev[TauType::tau23]->array(mfi);
506  auto t32_arr = (Tau_lev[TauType::tau32]) ? Tau_lev[TauType::tau32]->array(mfi) : Array4<Real>{};
507 
508  auto hfx3_arr = zheat_flux->array(mfi);
509  auto qfx3_arr = (zqv_flux) ? zqv_flux->array(mfi) : Array4<Real>{};
510 
511  // Rotated stress vars
512  auto t11_arr = (m_rotate) ? Tau_lev[TauType::tau11]->array(mfi) : Array4<Real>{};
513  auto t22_arr = (m_rotate) ? Tau_lev[TauType::tau22]->array(mfi) : Array4<Real>{};
514  auto t33_arr = (m_rotate) ? Tau_lev[TauType::tau33]->array(mfi) : Array4<Real>{};
515  auto t12_arr = (m_rotate) ? Tau_lev[TauType::tau12]->array(mfi) : Array4<Real>{};
516  auto t21_arr = (m_rotate) ? Tau_lev[TauType::tau21]->array(mfi) : Array4<Real>{};
517 
518  auto hfx1_arr = (m_rotate) ? xheat_flux->array(mfi) : Array4<Real>{};
519  auto hfx2_arr = (m_rotate) ? yheat_flux->array(mfi) : Array4<Real>{};
520  auto qfx1_arr = (m_rotate && xqv_flux) ? xqv_flux->array(mfi) : Array4<Real>{};
521  auto qfx2_arr = (m_rotate && yqv_flux) ? yqv_flux->array(mfi) : Array4<Real>{};
522 
523  // Terrain
524  const auto zphys_arr = (z_phys) ? z_phys->const_array(mfi) : Array4<const Real>{};
525 
526  // Get average arrays
527  const auto *const u_mean = m_ma.get_average(lev,0);
528  const auto *const v_mean = m_ma.get_average(lev,1);
529  const auto *const t_mean = m_ma.get_average(lev,2);
530  const auto *const q_mean = m_ma.get_average(lev,3);
531  const auto *const u_mag_mean = m_ma.get_average(lev,5);
532 
533  const auto um_arr = u_mean->array(mfi);
534  const auto vm_arr = v_mean->array(mfi);
535  const auto tm_arr = t_mean->array(mfi);
536  const auto qm_arr = q_mean->array(mfi);
537  const auto umm_arr = u_mag_mean->array(mfi);
538 
539  // Get derived arrays
540  const auto u_star_arr = u_star[lev]->array(mfi);
541  const auto t_star_arr = t_star[lev]->array(mfi);
542  const auto q_star_arr = q_star[lev]->array(mfi);
543  const auto t_surf_arr = t_surf[lev]->array(mfi);
544  const auto q_surf_arr = q_surf[lev]->array(mfi);
545 
546  // Get LSM fluxes
547  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
548  Array4<int> {};
549  auto lsm_t_flux_arr = Array4<Real> {};
550  auto lsm_q_flux_arr = Array4<Real> {};
551  auto lsm_tau13_arr = Array4<Real> {};
552  auto lsm_tau23_arr = Array4<Real> {};
553  for (int n(0); n<m_lsm_flux_lev[lev].size(); ++n) {
554  if (toLower(m_lsm_flux_name[n]) == "t_flux") { lsm_t_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
555  if (toLower(m_lsm_flux_name[n]) == "q_flux") { lsm_q_flux_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
556  if (toLower(m_lsm_flux_name[n]) == "tau13") { lsm_tau13_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
557  if (toLower(m_lsm_flux_name[n]) == "tau23") { lsm_tau23_arr = m_lsm_flux_lev[lev][n]->array(mfi); }
558  }
559 
560 
561  // Rho*Theta flux
562  //============================================================================
563  Box bx = mfi.tilebox();
564  if (bx.smallEnd(2) != klo) { continue; }
565  bx.makeSlab(2,klo);
566  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
567  {
568  // Valid theta flux from LSM and over land. The LSM writes the
569  // lsm_flux_undefined sentinel for cells it did not process (sea-ice /
570  // open water); fall back to MOST there instead of applying garbage.
571  Real Tflux;
572  int is_land = (lmask_arr) ? lmask_arr(i,j,0) : 1;
573  if (lsm_t_flux_arr && is_land && lsm_t_flux_arr(i,j,0) < lsm_flux_undefined) {
574  Tflux = lsm_t_flux_arr(i,j,0);
575  } else {
576  Tflux = flux_comp.compute_t_flux(i, j, k,
577  cons_arr, velx_arr, vely_arr,
578  umm_arr, tm_arr, u_star_arr,
579  t_star_arr, t_surf_arr);
580 
581  }
582 
583  // Do scalar flux rotations?
584  if (rotate) {
585  rotate_scalar_flux(i, j, k, Tflux, dxInv, zphys_arr,
586  hfx1_arr, hfx2_arr, hfx3_arr);
587  } else {
588  hfx3_arr(i,j,klo) = Tflux;
589  }
590  });
591 
592  // Rho*Qv flux
593  //============================================================================
594  if (use_moisture) {
595  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
596  {
597  // Valid qv flux from LSM and over land (sentinel -> fall back to MOST)
598  Real Qflux;
599  int is_land = (lmask_arr) ? lmask_arr(i,j,0) : 1;
600  if (lsm_q_flux_arr && is_land && lsm_q_flux_arr(i,j,0) < lsm_flux_undefined) {
601  Qflux = lsm_q_flux_arr(i,j,0);
602  } else {
603  Qflux = flux_comp.compute_q_flux(i, j, k,
604  cons_arr, velx_arr, vely_arr,
605  umm_arr, qm_arr, u_star_arr,
606  q_star_arr, q_surf_arr);
607  }
608 
609  // Do scalar flux rotations?
610  if (rotate) {
611  rotate_scalar_flux(i, j, k, Qflux, dxInv, zphys_arr,
612  qfx1_arr, qfx2_arr, qfx3_arr);
613  } else {
614  qfx3_arr(i,j,k) = Qflux;
615  }
616  });
617  } // custom
618 
619  if (!rotate) {
620  // Rho*u flux
621  //============================================================================
622  Box bxx = surroundingNodes(bx,0);
623  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
624  {
625  // Valid tau13 from LSM and over land. A side that is land but
626  // whose LSM flux is the sentinel (sea-ice / open water) is treated
627  // as non-LSM so that side uses the MOST stress instead.
628  Real stressx;
629  int is_land_hi = ((lmask_arr) ? lmask_arr(i ,j,0) : 1)
630  && (!lsm_tau13_arr || lsm_tau13_arr(i ,j,0) < lsm_flux_undefined);
631  int is_land_lo = ((lmask_arr) ? lmask_arr(i-1,j,0) : 1)
632  && (!lsm_tau13_arr || lsm_tau13_arr(i-1,j,0) < lsm_flux_undefined);
633  if (lsm_tau13_arr && (is_land_hi || is_land_lo)) {
634  stressx = zero;
635  if (!is_land_hi || !is_land_lo) {
636  stressx += myhalf * flux_comp.compute_u_flux(i, j, k,
637  cons_arr, velx_arr, vely_arr,
638  umm_arr, um_arr, u_star_arr);
639  }
640  if (is_land_hi) {
641  stressx += myhalf * lsm_tau13_arr(i ,j,0);
642  }
643  if (is_land_lo) {
644  stressx += myhalf * lsm_tau13_arr(i-1,j,0);
645  }
646  } else {
647  stressx = flux_comp.compute_u_flux(i, j, k,
648  cons_arr, velx_arr, vely_arr,
649  umm_arr, um_arr, u_star_arr);
650  }
651 
652  t13_arr(i,j,k) = stressx;
653  if (t31_arr) { t31_arr(i,j,k) = stressx; }
654  });
655 
656  // Rho*v flux
657  //============================================================================
658  Box bxy = surroundingNodes(bx,1);
659  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
660  {
661  // Valid tau23 from LSM and over land (sentinel side -> MOST stress)
662  Real stressy;
663  int is_land_hi = ((lmask_arr) ? lmask_arr(i,j ,0) : 1)
664  && (!lsm_tau23_arr || lsm_tau23_arr(i,j ,0) < lsm_flux_undefined);
665  int is_land_lo = ((lmask_arr) ? lmask_arr(i,j-1,0) : 1)
666  && (!lsm_tau23_arr || lsm_tau23_arr(i,j-1,0) < lsm_flux_undefined);
667  if (lsm_tau23_arr && (is_land_hi || is_land_lo)) {
668  stressy = zero;
669  if (!is_land_hi || !is_land_lo) {
670  stressy += myhalf * flux_comp.compute_v_flux(i, j, k,
671  cons_arr, velx_arr, vely_arr,
672  umm_arr, vm_arr, u_star_arr);
673  }
674  if (is_land_hi) {
675  stressy += myhalf * lsm_tau23_arr(i,j ,0);
676  }
677  if (is_land_lo) {
678  stressy += myhalf * lsm_tau23_arr(i,j-1,0);
679  }
680  } else {
681  stressy = flux_comp.compute_v_flux(i, j, k,
682  cons_arr, velx_arr, vely_arr,
683  umm_arr, vm_arr, u_star_arr);
684  }
685 
686  t23_arr(i,j,k) = stressy;
687  if (t32_arr) { t32_arr(i,j,k) = stressy; }
688  });
689  } else {
690  // All fluxes with rotation
691  //============================================================================
692  Box bxxy = convert(bx, IntVect(1,1,0));
693  ParallelFor(bxxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
694  {
695  Real stresst = flux_comp.compute_u_flux(i, j, k,
696  cons_arr, velx_arr, vely_arr,
697  umm_arr, um_arr, u_star_arr);
698  rotate_stress_tensor(i, j, k, stresst, dxInv, zphys_arr,
699  velx_arr, vely_arr, velz_arr,
700  t11_arr, t22_arr, t33_arr,
701  t12_arr, t21_arr,
702  t13_arr, t31_arr,
703  t23_arr, t32_arr);
704  });
705  }
706  } // mfiter
707 }
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau32
Definition: ERF_DataStruct.H:32
@ tau31
Definition: ERF_DataStruct.H:32
@ tau21
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
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:611
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:632
@ xvel
Definition: ERF_IndexDefines.H:175
@ zvel
Definition: ERF_IndexDefines.H:177
@ yvel
Definition: ERF_IndexDefines.H:176
Here is the call graph for this function:

◆ compute_SurfaceLayer_bcs_EB() [1/2]

template<typename FluxCalc >
void SurfaceLayer::compute_SurfaceLayer_bcs_EB ( const int &  lev,
amrex::Vector< const amrex::MultiFab * >  mfs,
amrex::Vector< 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 FluxCalc &  flux_comp 
)

◆ compute_SurfaceLayer_bcs_EB() [2/2]

template<typename FluxCalc >
void SurfaceLayer::compute_SurfaceLayer_bcs_EB ( const int &  lev,
Vector< const MultiFab * >  mfs,
Vector< Vector< std::unique_ptr< MultiFab >>> &  Tau_EB,
[[maybe_unused] ] MultiFab *  xheat_flux,
[[maybe_unused] ] MultiFab *  yheat_flux,
MultiFab *  Hfx3_EB,
[[maybe_unused] ] MultiFab *  xqv_flux,
[[maybe_unused] ] MultiFab *  yqv_flux,
[[maybe_unused] ] MultiFab *  zqv_flux,
const FluxCalc &  flux_comp 
)

Function to calculate MOST fluxes for EB.

Parameters
[in]levCurrent level
[in,out]mfsMultiFabs to populate
[in]eddyDiffsDiffusion coefficients from turbulence model
[in]flux_compstructure to compute fluxes
729 {
730  // Get EB flags for all centerings
731  const auto& cc_flags = m_eb_vec[lev]->get_const_factory()->getMultiEBCellFlagFab();
732  const auto& u_flags = m_eb_vec[lev]->get_u_const_factory()->getMultiEBCellFlagFab();
733  const auto& v_flags = m_eb_vec[lev]->get_v_const_factory()->getMultiEBCellFlagFab();
734  const auto& w_flags = m_eb_vec[lev]->get_w_const_factory()->getMultiEBCellFlagFab();
735 
736  const auto& cc_vfrac = m_eb_vec[lev]->get_const_factory()->getVolFrac();
737  const auto& u_vfrac = m_eb_vec[lev]->get_u_const_factory()->getVolFrac();
738  const auto& v_vfrac = m_eb_vec[lev]->get_v_const_factory()->getVolFrac();
739  const auto& w_vfrac = m_eb_vec[lev]->get_w_const_factory()->getVolFrac();
740  const auto& cc_bnorm = m_eb_vec[lev]->get_const_factory()->getBndryNormal();
741  const auto& u_bnorm = m_eb_vec[lev]->get_u_const_factory()->getBndryNorm();
742  const auto& v_bnorm = m_eb_vec[lev]->get_v_const_factory()->getBndryNorm();
743  const auto& w_bnorm = m_eb_vec[lev]->get_w_const_factory()->getBndryNorm();
744 
745  for (MFIter mfi(*mfs[0]); mfi.isValid(); ++mfi)
746  {
747  // Get flags for this box (all centerings)
748  const auto& cc_flag = cc_flags[mfi];
749  const auto& u_flag = u_flags[mfi];
750  const auto& v_flag = v_flags[mfi];
751  const auto& w_flag = w_flags[mfi];
752 
753  // Skip boxes that have no cut cells at any centering
754  if (cc_flag.getType() != FabType::singlevalued &&
755  u_flag.getType() != FabType::singlevalued &&
756  v_flag.getType() != FabType::singlevalued &&
757  w_flag.getType() != FabType::singlevalued
758  ) continue;
759 
760  // Get EB flag and volfrac arrays
761  auto const cc_flag_arr = cc_flag.const_array();
762  auto const u_flag_arr = u_flag.const_array();
763  auto const v_flag_arr = v_flag.const_array();
764  auto const w_flag_arr = w_flag.const_array();
765 
766  auto const cc_vfrac_arr = cc_vfrac.const_array(mfi);
767  auto const u_vfrac_arr = u_vfrac.const_array(mfi);
768  auto const v_vfrac_arr = v_vfrac.const_array(mfi);
769  auto const w_vfrac_arr = w_vfrac.const_array(mfi);
770  auto const bnorm_arr = cc_bnorm.const_array(mfi);
771  auto const u_bnorm_arr = u_bnorm.const_array(mfi);
772  auto const v_bnorm_arr = v_bnorm.const_array(mfi);
773  auto const w_bnorm_arr = w_bnorm.const_array(mfi);
774 
775  // Get field arrays
776  const auto cons_arr = mfs[Vars::cons]->array(mfi);
777  const auto velx_arr = mfs[Vars::xvel]->array(mfi);
778  const auto vely_arr = mfs[Vars::yvel]->array(mfi);
779  const auto velz_arr = mfs[Vars::zvel]->array(mfi);
780 
781  // Diffusive stress vars - t13 and t23 components for all grid types
782  auto u_t13_arr = Tau_EB[EBTauType::tau_eb13][EBGridType::xface]->array(mfi);
783  auto v_t13_arr = Tau_EB[EBTauType::tau_eb13][EBGridType::yface]->array(mfi);
784  auto w_t13_arr = Tau_EB[EBTauType::tau_eb13][EBGridType::zface]->array(mfi);
785 
786  auto u_t23_arr = Tau_EB[EBTauType::tau_eb23][EBGridType::xface]->array(mfi);
787  auto v_t23_arr = Tau_EB[EBTauType::tau_eb23][EBGridType::yface]->array(mfi);
788  auto w_t23_arr = Tau_EB[EBTauType::tau_eb23][EBGridType::zface]->array(mfi);
789 
790  auto hfx3_arr = Hfx3_EB->array(mfi);
791 
792  // Get average arrays
793  const auto *const u_mean = m_ma.get_average(lev,0);
794  const auto *const v_mean = m_ma.get_average(lev,1);
795  const auto *const t_mean = m_ma.get_average(lev,2);
796  // const auto *const q_mean = m_ma.get_average(lev,3);
797  const auto *const u_mag_mean = m_ma.get_average(lev,5);
798 
799  const auto um_arr = u_mean->array(mfi);
800  const auto vm_arr = v_mean->array(mfi);
801  const auto tm_arr = t_mean->array(mfi);
802  // const auto qm_arr = q_mean->array(mfi);
803  const auto umm_arr = u_mag_mean->array(mfi);
804 
805  // Get derived arrays
806  const auto u_star_arr = u_star[lev]->array(mfi);
807  const auto t_star_arr = t_star[lev]->array(mfi);
808  const auto t_surf_arr = t_surf[lev]->array(mfi);
809 
810  // Rho*Theta flux
811  //============================================================================
812  Box bx = mfi.tilebox();
813  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
814  {
815  if (cc_flag_arr(i,j,k).isSingleValued()) {
816  Real Tflux = flux_comp.compute_t_flux(i, j, k,
817  cons_arr, velx_arr, vely_arr, velz_arr,
818  umm_arr, tm_arr, u_star_arr,
819  t_star_arr, t_surf_arr,
820  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
821  bnorm_arr);
822  hfx3_arr(i,j,k) = Tflux;
823  }
824  });
825 
826  // Rho*u flux
827  //============================================================================
828  Box bxx = surroundingNodes(bx,0);
829  Box bxy = surroundingNodes(bx,1);
830  Box bxz = surroundingNodes(bx,2);
831  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
832  {
833  if (u_flag_arr(i,j,k).isSingleValued()) {
834  Real stressx = flux_comp.compute_u_flux(i, j, k,
835  cons_arr, velx_arr, vely_arr, velz_arr,
836  umm_arr, um_arr, u_star_arr,
837  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
838  cc_vfrac_arr, cc_flag_arr,
839  u_bnorm_arr, 0);
840  u_t13_arr(i,j,k) = stressx;
841  }
842  });
843  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
844  {
845  if (v_flag_arr(i,j,k).isSingleValued()) {
846  Real stressx = flux_comp.compute_u_flux(i, j, k,
847  cons_arr, velx_arr, vely_arr, velz_arr,
848  umm_arr, um_arr, u_star_arr,
849  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
850  cc_vfrac_arr, cc_flag_arr,
851  v_bnorm_arr, 1);
852  v_t13_arr(i,j,k) = stressx;
853  }
854  });
855  ParallelFor(bxz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
856  {
857  if (u_flag_arr(i,j,k).isSingleValued()) {
858  Real stressx = flux_comp.compute_u_flux(i, j, k,
859  cons_arr, velx_arr, vely_arr, velz_arr,
860  umm_arr, um_arr, u_star_arr,
861  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
862  cc_vfrac_arr, cc_flag_arr,
863  w_bnorm_arr, 2);
864  w_t13_arr(i,j,k) = stressx;
865  }
866  });
867 
868  // Rho*v flux
869  //============================================================================
870  ParallelFor(bxx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
871  {
872  if (u_flag_arr(i,j,k).isSingleValued()) {
873  Real stressy = flux_comp.compute_v_flux(i, j, k,
874  cons_arr, velx_arr, vely_arr, velz_arr,
875  umm_arr, vm_arr, u_star_arr,
876  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
877  cc_vfrac_arr, cc_flag_arr, u_bnorm_arr, 0);
878  u_t23_arr(i,j,k) = stressy;
879  }
880  });
881  ParallelFor(bxy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
882  {
883  if (v_flag_arr(i,j,k).isSingleValued()) {
884  Real stressy = flux_comp.compute_v_flux(i, j, k,
885  cons_arr, velx_arr, vely_arr, velz_arr,
886  umm_arr, vm_arr, u_star_arr,
887  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
888  cc_vfrac_arr, cc_flag_arr, v_bnorm_arr, 1);
889  v_t23_arr(i,j,k) = stressy;
890  }
891  });
892  ParallelFor(bxz, [=] AMREX_GPU_DEVICE (int i, int j, int k)
893  {
894  if (w_flag_arr(i,j,k).isSingleValued()) {
895  Real stressy = flux_comp.compute_v_flux(i, j, k,
896  cons_arr, velx_arr, vely_arr, velz_arr,
897  umm_arr, vm_arr, u_star_arr,
898  u_vfrac_arr, v_vfrac_arr, w_vfrac_arr,
899  cc_vfrac_arr, cc_flag_arr, w_bnorm_arr, 2);
900  w_t23_arr(i,j,k) = stressy;
901  }
902  });
903  } // mfiter
904 }
@ tau_eb23
Definition: ERF_EBStruct.H:16
@ tau_eb13
Definition: ERF_EBStruct.H:16
@ yface
Definition: ERF_EBStruct.H:20
@ zface
Definition: ERF_EBStruct.H:20
@ xface
Definition: ERF_EBStruct.H:20
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 
)
1066 {
1067  // NOTE: We have already tested a moisture model exists
1068 
1069  // Populate q_surf with qsat over water
1070  auto dz = m_geom[lev].CellSize(2);
1071  const int klo = m_geom[lev].Domain().smallEnd(2);
1072  for (MFIter mfi(*q_surf[lev]); mfi.isValid(); ++mfi)
1073  {
1074  Box gtbx = mfi.growntilebox();
1075 
1076  if (gtbx.smallEnd(2) != klo) { continue; }
1077 
1078  auto t_surf_arr = t_surf[lev]->array(mfi);
1079  auto q_surf_arr = q_surf[lev]->array(mfi);
1080  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
1081  Array4<int> {};
1082  const auto cons_arr = cons_in.const_array(mfi);
1083  const auto z_arr = (z_phys_nd) ? z_phys_nd->const_array(mfi) :
1084  Array4<const Real> {};
1085 
1086  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1087  {
1088  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
1089  if (!is_land) {
1090  auto deltaZ = (z_arr) ? Compute_Zrel_AtCellCenter(i,j,k,z_arr) :
1091  myhalf*dz;
1092  auto Rho = cons_arr(i,j,k,Rho_comp);
1093  auto RTh = cons_arr(i,j,k,RhoTheta_comp);
1094  auto Qv = cons_arr(i,j,k,RhoQ1_comp) / Rho;
1095  auto P_cc = getPgivenRTh(RTh, Qv);
1096  P_cc += Rho*CONST_GRAV*deltaZ;
1097  P_cc *= Real(0.01);
1098  erf_qsatw(t_surf_arr(i,j,k), P_cc, q_surf_arr(i,j,k));
1099  }
1100  });
1101  }
1102  q_surf[lev]->FillBoundary(m_geom[lev].periodicity());
1103 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void erf_qsatw(amrex::Real t, amrex::Real p, amrex::Real &qsatw)
Definition: ERF_MicrophysicsUtils.H:228
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:389
@ dz
Definition: ERF_AdvanceWSM6.cpp:104
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 
)
978 {
979  int n_times_in_sst = static_cast<int>(m_sst_lev[lev].size());
980 
982 
983  int n_time_lo, n_time_hi;
984  Real alpha;
985 
986  if (n_times_in_sst > 1) {
987  n_time_lo = static_cast<int>( elapsed_time_since_start_low / dT);
988  alpha = (elapsed_time_since_start_low - n_time_lo * dT) / dT;
989 
991 
992  n_time_hi = n_time_lo + 1;
993 
994  // Do not over run the last sst file
995  if (m_start_low_time + elapsed_time_since_start_low >= m_final_low_time) {
996  n_time_lo = static_cast<int>(m_sst_lev[lev].size())-1;
997  n_time_hi = n_time_lo;
998  alpha = zero;
999  }
1000 
1001  AMREX_ALWAYS_ASSERT( (n_time_lo >= 0) && (n_time_hi < m_sst_lev[lev].size()));
1002  } else {
1003  n_time_lo = 0;
1004  n_time_hi = 0;
1005  alpha = one;
1006  }
1007  AMREX_ALWAYS_ASSERT( alpha >= zero && alpha <= one);
1008 
1009  Real oma = one - alpha;
1010 
1011  // Define a default land surface temperature if we don't read in tsk
1013 
1014  bool use_tsk = (m_tsk_lev[lev][0]);
1015  bool ignore_sst = m_ignore_sst;
1016 
1017  const int klo = m_geom[lev].Domain().smallEnd(2);
1018 
1019  // Populate t_surf
1020  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
1021  {
1022  Box gtbx = mfi.growntilebox();
1023 
1024  if (gtbx.smallEnd(2) != klo) { continue; }
1025 
1026  auto t_surf_arr = t_surf[lev]->array(mfi);
1027 
1028  const auto sst_lo_arr = m_sst_lev[lev][n_time_lo]->const_array(mfi);
1029  const auto sst_hi_arr = m_sst_lev[lev][n_time_hi]->const_array(mfi);
1030 
1031  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
1032  Array4<int> {};
1033 
1034  if (use_tsk) {
1035  const auto tsk_arr = m_tsk_lev[lev][n_time_lo]->const_array(mfi);
1036  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1037  {
1038  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
1039  if (!is_land && !ignore_sst) {
1040  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
1041  + alpha * sst_hi_arr(i,j,k);
1042  } else {
1043  t_surf_arr(i,j,k) = tsk_arr(i,j,k);
1044  }
1045  });
1046  } else {
1047  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1048  {
1049  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
1050  if (!is_land) {
1051  t_surf_arr(i,j,k) = oma * sst_lo_arr(i,j,k)
1052  + alpha * sst_hi_arr(i,j,k);
1053  } else {
1054  t_surf_arr(i,j,k) = lst;
1055  }
1056  });
1057  }
1058  }
1059  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
1060 }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_sst_lev
Definition: ERF_SurfaceLayer.H:771
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_tsk_lev
Definition: ERF_SurfaceLayer.H:772
real(kind=kind_phys), parameter, public alpha
Definition: ERF_module_mp_wsm6.F90:44
Here is the call graph for this function:

◆ get_lmask()

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

◆ get_lsm_tsurf()

void SurfaceLayer::get_lsm_tsurf ( const int &  lev)
1107 {
1108  const int klo = m_geom[lev].Domain().smallEnd(2);
1109  const bool has_sea_tsurf = (m_has_ocean_lsm_tsurf &&
1110  amrex::toLower(m_lsm_data_name[m_lsm_tsurf_indx]) == "t_surf");
1111  for (MFIter mfi(*t_surf[lev]); mfi.isValid(); ++mfi)
1112  {
1113  Box gtbx = mfi.growntilebox();
1114 
1115  if (gtbx.smallEnd(2) != klo) { continue; }
1116 
1117  // NOTE: LSM does not carry lateral ghost cells.
1118  // This copies the valid box into the ghost cells.
1119  // Fillboundary is called after this to pick up the
1120  // interior ghost and periodic directions.
1121  Box vbx = mfi.validbox();
1122  int i_lo = vbx.smallEnd(0); int i_hi = vbx.bigEnd(0);
1123  int j_lo = vbx.smallEnd(1); int j_hi = vbx.bigEnd(1);
1124 
1125  auto t_surf_arr = t_surf[lev]->array(mfi);
1126  auto lmask_arr = (m_lmask_lev[lev][0]) ? m_lmask_lev[lev][0]->array(mfi) :
1127  Array4<int> {};
1128  const auto lsm_arr = m_lsm_data_lev[lev][m_lsm_tsurf_indx]->const_array(mfi);
1129 
1130  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1131  {
1132  int is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
1133  if ((!has_sea_tsurf && is_land) ||
1134  (has_sea_tsurf && !is_land)) {
1135  int li = amrex::min(amrex::max(i, i_lo), i_hi);
1136  int lj = amrex::min(amrex::max(j, j_lo), j_hi);
1137  t_surf_arr(i,j,k) = lsm_arr(li,lj,k);
1138  }
1139  });
1140  }
1141 }
amrex::Vector< std::string > m_lsm_data_name
Definition: ERF_SurfaceLayer.H:776
int m_lsm_tsurf_indx
Definition: ERF_SurfaceLayer.H:751
amrex::Vector< amrex::Vector< amrex::MultiFab * > > m_lsm_data_lev
Definition: ERF_SurfaceLayer.H:774
bool m_has_ocean_lsm_tsurf
Definition: ERF_SurfaceLayer.H:750
Here is the call graph for this function:

◆ get_mac_avg()

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

◆ get_olen()

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

◆ get_pblh()

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

◆ get_q_star()

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

◆ get_q_surf()

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

◆ get_t_star()

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

◆ get_t_surf()

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

◆ get_u_star()

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

◆ get_w_star()

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

◆ get_z0()

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

◆ get_zref()

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

◆ have_variable_sea_roughness()

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

◆ 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
397 {
399  moeng_flux flux_comp;
400  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
401  xheat_flux, yheat_flux, zheat_flux,
402  xqv_flux, yqv_flux, zqv_flux,
403  z_phys, flux_comp);
404  } else if (flux_type == FluxCalcType::DONELAN) {
405  donelan_flux flux_comp;
406  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
407  xheat_flux, yheat_flux, zheat_flux,
408  xqv_flux, yqv_flux, zqv_flux,
409  z_phys, flux_comp);
410  } else if (flux_type == FluxCalcType::ROTATE) {
411  rotate_flux flux_comp;
412  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
413  xheat_flux, yheat_flux, zheat_flux,
414  xqv_flux, yqv_flux, zqv_flux,
415  z_phys, flux_comp);
416  } else if (flux_type == FluxCalcType::RICO) {
418  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
419  xheat_flux, yheat_flux, zheat_flux,
420  xqv_flux, yqv_flux, zqv_flux,
421  z_phys, flux_comp);
422  } else if (flux_type == FluxCalcType::BULK_COEFF) {
423  bulk_coeff_flux flux_comp(m_Cd, m_Ch, m_Cq);
424  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
425  xheat_flux, yheat_flux, zheat_flux,
426  xqv_flux, yqv_flux, zqv_flux,
427  z_phys, flux_comp);
428  } else if (flux_type == FluxCalcType::CUSTOM) {
429  custom_flux flux_comp(specified_rho_surf);
430  compute_SurfaceLayer_bcs(lev, mfs, Tau_lev,
431  xheat_flux, yheat_flux, zheat_flux,
432  xqv_flux, yqv_flux, zqv_flux,
433  z_phys, flux_comp);
434  } else {
435  amrex::Abort("Unknown surface layer flux calculation type");
436  }
437 }
bool specified_rho_surf
Definition: ERF_SurfaceLayer.H:737
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:2217
Definition: ERF_MOSTStress.H:2096
Definition: ERF_MOSTStress.H:1963
Definition: ERF_MOSTStress.H:1796
Definition: ERF_MOSTStress.H:2344
Definition: ERF_MOSTStress.H:2476

◆ impose_SurfaceLayer_bcs_EB()

void SurfaceLayer::impose_SurfaceLayer_bcs_EB ( const int &  lev,
amrex::Vector< const amrex::MultiFab * >  mfs,
amrex::Vector< 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 
)

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
456 {
458  moeng_flux_eb flux_comp;
459  compute_SurfaceLayer_bcs_EB(lev, mfs, Tau_EB,
460  xheat_flux, yheat_flux, Hfx3_EB,
461  xqv_flux, yqv_flux, zqv_flux,
462  flux_comp);
463  } else {
464  amrex::Abort("Not implemented surface layer flux calculation type for EB");
465  }
466 }
void compute_SurfaceLayer_bcs_EB(const int &lev, amrex::Vector< const amrex::MultiFab * > mfs, amrex::Vector< 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 FluxCalc &flux_comp)
Definition: ERF_EBMOSTStress.H:298

◆ init_tke_from_ustar()

void SurfaceLayer::init_tke_from_ustar ( const int &  lev,
amrex::MultiFab &  cons,
const std::unique_ptr< amrex::MultiFab > &  z_phys_nd,
const amrex::Real  tkefac = one,
const amrex::Real  zscale = amrex::Real(700.0) 
)
1176 {
1177  Print() << "Initializing TKE from surface layer ustar on level " << lev << std::endl;
1178 
1179  // Handle vertical decomposition by selectively copying into
1180  // a FArrayBox section on each rank. Then doing a reduce real sum
1181  // and broadcasting to each rank. No mask since all CC data
1182  const int klo = m_geom[lev].Domain().smallEnd(2);
1183  Box bx_lo = u_star[lev]->boxArray().minimalBox();
1184  FArrayBox u_star_lo(bx_lo, 1); u_star_lo.setVal<RunOn::Device>(0);
1185  FArrayBox z_surf_lo(bx_lo, 1); z_surf_lo.setVal<RunOn::Device>(0);
1186  Real* ustar_ptr = u_star_lo.dataPtr();
1187  Real* zsurf_ptr = z_surf_lo.dataPtr();
1188  for (MFIter mfi(cons); mfi.isValid(); ++mfi)
1189  {
1190  Box vbx = mfi.validbox();
1191  if (vbx.smallEnd(2) != klo) { continue; }
1192  vbx.makeSlab(2,0);
1193 
1194  auto const& u_star_arr = u_star[lev]->const_array(mfi);
1195  auto u_star_all = u_star_lo.array();
1196 
1197  auto const& z_phys_arr = z_phys_nd->const_array(mfi);
1198  auto z_surf_all = z_surf_lo.array();
1199 
1200  ParallelFor(vbx, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
1201  {
1202  u_star_all(i,j,0) = u_star_arr(i,j,0);
1203  z_surf_all(i,j,0) = fourth * ( z_phys_arr(i ,j ,klo) + z_phys_arr(i+1,j ,klo)
1204  + z_phys_arr(i ,j+1,klo) + z_phys_arr(i+1,j+1,klo) );
1205  });
1206  }
1207  ParallelDescriptor::ReduceRealSum(ustar_ptr, static_cast<int>(bx_lo.numPts()));
1208  ParallelDescriptor::ReduceRealSum(zsurf_ptr, static_cast<int>(bx_lo.numPts()));
1209 
1210  // Now work on all boxes (ustar has been filled above)
1211  constexpr Real small = Real(0.01);
1212  for (MFIter mfi(cons); mfi.isValid(); ++mfi)
1213  {
1214  Box vbx = mfi.validbox();
1215 
1216  auto const& u_star_arr = u_star_lo.const_array();
1217  auto const& z_surf_arr = z_surf_lo.const_array();
1218  auto const& z_phys_arr = z_phys_nd->const_array(mfi);
1219 
1220  auto const& cons_arr = cons.array(mfi);
1221 
1222  ParallelFor(vbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
1223  {
1224  Real rho = cons_arr(i, j, k, Rho_comp);
1225  Real ust = u_star_arr(i, j, 0);
1226  Real tke0 = tkefac * ust * ust; // surface value
1227  Real zagl = Compute_Z_AtCellCenter(i, j, k, z_phys_arr) - z_surf_arr(i,j,0);
1228 
1229  // linearly tapering profile -- following WRF, approximate top of
1230  // PBL as ustar * zscale
1231  cons_arr(i, j, k, RhoKE_comp) = rho * tke0 * std::max(
1232  (ust * zscale - zagl) / (std::max(ust, small) * zscale),
1233  small);
1234  });
1235  }
1236 }
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtCellCenter(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:361
Here is the call graph for this function:

◆ lmask_min_reduce()

int SurfaceLayer::lmask_min_reduce ( amrex::iMultiFab &  lmask,
const int &  nghost 
)
inline
652  {
653  int lmask_min = amrex::ReduceMin(lmask, nghost, [=] AMREX_GPU_HOST_DEVICE(
654  amrex::Box const& bx, amrex::Array4<int const> const& lm_arr) -> int
655  {
656  int locmin = std::numeric_limits<int>::max();
657  const auto lo = lbound(bx);
658  const auto hi = ubound(bx);
659  for (int j = lo.y; j <= hi.y; ++j) {
660  for (int i = lo.x; i <= hi.x; ++i) {
661  locmin = std::min(locmin, lm_arr(i, j, 0));
662  }
663  }
664  return locmin;
665  });
666 
667  return lmask_min;
668  }

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

◆ read_custom_roughness()

void SurfaceLayer::read_custom_roughness ( const int &  lev,
const std::string &  fname 
)
1242 {
1243  // Read the file if we have it
1244  if (!fname.empty()) {
1245  // Only the ioproc reads the file
1246  Gpu::HostVector<Real> m_x,m_y,m_z0;
1247  if (ParallelDescriptor::IOProcessor()) {
1248  Print()<<"Reading MOST roughness file at level " << lev << " : " << fname << std::endl;
1249  std::ifstream file(fname);
1250  Real value1,value2,value3;
1251  while(file>>value1>>value2>>value3){
1252  m_x.push_back(value1);
1253  m_y.push_back(value2);
1254  m_z0.push_back(value3);
1255  }
1256  file.close();
1257 
1258  AMREX_ALWAYS_ASSERT(m_x.size() == m_y.size());
1259  AMREX_ALWAYS_ASSERT(m_x.size() == m_z0.size());
1260  }
1261 
1262  // Broadcast the whole domain to every rank
1263  int ioproc = ParallelDescriptor::IOProcessorNumber();
1264  int nnode = static_cast<int>(m_x.size());
1265  ParallelDescriptor::Bcast(&nnode, 1, ioproc);
1266 
1267  if (!ParallelDescriptor::IOProcessor()) {
1268  m_x.resize(nnode);
1269  m_y.resize(nnode);
1270  m_z0.resize(nnode);
1271  }
1272  ParallelDescriptor::Bcast(m_x.data() , nnode, ioproc);
1273  ParallelDescriptor::Bcast(m_y.data() , nnode, ioproc);
1274  ParallelDescriptor::Bcast(m_z0.data(), nnode, ioproc);
1275 
1276  // Copy data to the GPU
1277  Gpu::DeviceVector<Real> d_x(nnode),d_y(nnode),d_z0(nnode);
1278  Gpu::copy(Gpu::hostToDevice, m_x.begin(), m_x.end(), d_x.begin());
1279  Gpu::copy(Gpu::hostToDevice, m_y.begin(), m_y.end(), d_y.begin());
1280  Gpu::copy(Gpu::hostToDevice, m_z0.begin(), m_z0.end(), d_z0.begin());
1281  Real* xp = d_x.data();
1282  Real* yp = d_y.data();
1283  Real* z0p = d_z0.data();
1284 
1285  // Each rank populates it's z_0[lev] MultiFab
1286  const int klo = m_geom[lev].Domain().smallEnd(2);
1287  for (MFIter mfi(z_0[lev]); mfi.isValid(); ++mfi)
1288  {
1289  Box gtbx = mfi.growntilebox();
1290 
1291  if (gtbx.smallEnd(2) != klo) { continue; }
1292 
1293  // Populate z_phys data
1294  Real tol = Real(1.0e-4);
1295  auto dx = m_geom[lev].CellSizeArray();
1296  auto ProbLoArr = m_geom[lev].ProbLoArray();
1297  int ilo = m_geom[lev].Domain().smallEnd(0);
1298  int jlo = m_geom[lev].Domain().smallEnd(1);
1299  int ihi = m_geom[lev].Domain().bigEnd(0);
1300  int jhi = m_geom[lev].Domain().bigEnd(1);
1301 
1302  Array4<Real> const& z0_arr = z_0[lev].array(mfi);
1303  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE (int i, int j, int /*k*/)
1304  {
1305  // Clip indices for ghost-cells
1306  int ii = amrex::min(amrex::max(i,ilo),ihi);
1307  int jj = amrex::min(amrex::max(j,jlo),jhi);
1308 
1309  // Location of nodes
1310  Real x = ProbLoArr[0] + ii * dx[0];
1311  Real y = ProbLoArr[1] + jj * dx[1];
1312  int inode = ii + jj * (ihi-ilo+2); // stride is Nx+1
1313  if (std::sqrt(amrex::Math::powi<2>(x-xp[inode])+amrex::Math::powi<2>(y-yp[inode])) < tol) {
1314  z0_arr(i,j,klo) = z0p[inode];
1315  } else {
1316  // Unexpected list order, do brute force search
1317  Real z0loc = zero;
1318  bool found = false;
1319  for (int n=0; n<nnode; ++n) {
1320  Real delta=std::sqrt(amrex::Math::powi<2>(x-xp[n])+amrex::Math::powi<2>(y-yp[n]));
1321  if (delta < tol) {
1322  found = true;
1323  z0loc = z0p[n];
1324  break;
1325  }
1326  }
1327  AMREX_ASSERT_WITH_MESSAGE(found, "Location read from terrain file does not match the grid!");
1328  amrex::ignore_unused(found);
1329  z0_arr(i,j,klo) = z0loc;
1330  }
1331  });
1332  } // mfi
1333  } else {
1334  AMREX_ALWAYS_ASSERT(lev > 0);
1335 
1336  Print()<<"Interpolating MOST roughness at level " << lev << std::endl;
1337 
1338  // Create a BC mapper that uses FOEXTRAP at domain bndry
1339  Vector<int> bc_lo(3,ERFBCType::foextrap);
1340  Vector<int> bc_hi(3,ERFBCType::foextrap);
1341  Vector<BCRec> bcr; bcr.push_back(BCRec(bc_lo.data(),bc_hi.data()));
1342 
1343  // Create ref ratio
1344  IntVect ratio;
1345  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
1346  ratio[idim] = m_geom[lev].Domain().length(idim) / m_geom[0].Domain().length(idim);
1347  }
1348 
1349  // Create interp object and interpolate from the coarsest grid
1350  MFInterpolater* interp = &mf_cell_cons_interp;
1351  interp->interp(z_0[0] , 0,
1352  z_0[lev], 0,
1353  1, z_0[lev].nGrowVect(),
1354  m_geom[0], m_geom[lev],
1355  m_geom[lev].Domain(),ratio,
1356  bcr, 0);
1357  }
1358 }
@ m_y
Definition: ERF_DataStruct.H:25
@ m_x
Definition: ERF_DataStruct.H:24
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
@ foextrap
Definition: ERF_IndexDefines.H:242

Referenced by make_SurfaceLayer_at_level().

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

◆ set_q_surf()

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

◆ set_t_surf()

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

◆ update_fluxes()

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

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

Parameters
[in]levCurrent level
[in]max_itersmaximum iterations to use
19 {
20  // Update with SST/TSK data if we have a valid pointer
21  if (!m_has_ocean_lsm_tsurf &&
22  !m_sst_lev[lev].empty() && m_sst_lev[lev][0]) {
23  fill_tsurf_with_sst_and_tsk(lev, elapsed_time_since_start_low);
24  }
25 
26  // Apply heating rate if needed
29  update_surf_temp(elapsed_time_since_start_low);
30  }
31 
32  // Update qsurf with qsat over sea
33  if (use_moisture) {
34  fill_qsurf_with_qsat(lev, cons_in, z_phys_nd);
35  }
36 
37  // Update land surface temp if we have a valid pointer
38  if (m_has_lsm_tsurf) { get_lsm_tsurf(lev); }
39 
40  // Fill interior ghost cells
41  t_surf[lev]->FillBoundary(m_geom[lev].periodicity());
42 
43  // Compute plane averages for all vars (regardless of flux type)
45 
46  // NOTE: Do iterations to seed variables on the first step (LSM called post step)
47  //*******************************************************************************
48  // Update u*/T*/q*/L over land (iterations or from LSM fluxes)
49  if (m_has_lsm_fluxes && elapsed_time > zero) {
51  } else {
52  // ***************************************************************
53  // Iterate the fluxes if moeng type
54  // First iterate over land -- the only model for surface roughness
55  // over land is RoughCalcType::CONSTANT
56  // ***************************************************************
59  bool is_land = true;
60  // Do we have a constant flux for moisture over land?
61  bool cons_qflux = ( (moist_type == MoistCalcType::MOISTURE_FLUX) ||
63  if (m_terrain_type != TerrainType::EB) {
66  surface_flux most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
67  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
68  } else {
69  amrex::Abort("Unknown value for rough_type_land");
70  }
73  surface_temp most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
74  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
75  } else {
76  amrex::Abort("Unknown value for rough_type_land");
77  }
78  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
82  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
83  } else {
84  amrex::Abort("Unknown value for rough_type_land");
85  }
86  } else {
87  amrex::Abort("Unknown value for theta_type");
88  }
89  // EB
90  } else {
93  surface_flux_eb most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
94  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
95  } else {
96  amrex::Abort("Unknown value for rough_type_land");
97  }
100  surface_temp_eb most_flux(surf_temp_flux, surf_moist_flux, cons_qflux);
101  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
102  } else {
103  amrex::Abort("Unknown value for rough_type_land");
104  }
105  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
109  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
110  } else {
111  amrex::Abort("Unknown value for rough_type_land");
112  }
113  } else {
114  amrex::Abort("Unknown value for theta_type");
115  }
116  }
117  } // MOENG -- LAND
118  } // has_lsm_fluxes
119 
120  // ***************************************************************
121  // Iterate the fluxes if moeng type
122  // Next iterate over sea -- the models for surface roughness
123  // over sea are CHARNOCK, DONELAN, MODIFIED_CHARNOCK or WAVE_COUPLED
124  // NOTE: Sea surface fluxes are not supported for EB terrain
125  // ***************************************************************
126  if ((flux_type == FluxCalcType::MOENG ||
128  m_terrain_type != TerrainType::EB) {
129  bool is_land = false;
130  // NOTE: Do not allow default to adiabatic over sea (we have Qvs at surface)
131  // Do we have a constant flux for moisture over sea?
132  bool cons_qflux = (moist_type == MoistCalcType::MOISTURE_FLUX);
136  cnk_a, cnk_visc, cons_qflux);
137  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
140  depth, cons_qflux);
141  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
142  } else if (rough_type_sea == RoughCalcType::DONELAN) {
144  cons_qflux);
145  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
148  cons_qflux);
149  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
150  } else {
151  amrex::Abort("Unknown value for rough_type_sea");
152  }
153 
157  cnk_a, cnk_visc, cons_qflux);
158  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
161  depth, cons_qflux);
162  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
163  } else if (rough_type_sea == RoughCalcType::DONELAN) {
165  cons_qflux);
166  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
169  cons_qflux);
170  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
171  } else {
172  amrex::Abort("Unknown value for rough_type_sea");
173  }
174 
175  } else if ((theta_type == ThetaCalcType::ADIABATIC) &&
179  cnk_a, cnk_visc);
180  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
183  depth);
184  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
185  } else if (rough_type_sea == RoughCalcType::DONELAN) {
187  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
190  compute_fluxes(lev, max_iters, cons_in, most_flux, is_land);
191  } else {
192  amrex::Abort("Unknown value for rough_type_sea");
193  }
194  } else {
195  amrex::Abort("Unknown value for theta_type");
196  }
197  } // MOENG -- SEA
198 
200  if (custom_rhosurf > 0) {
201  specified_rho_surf = true;
202  u_star[lev]->setVal(std::sqrt(custom_rhosurf) * custom_ustar);
203  t_star[lev]->setVal(custom_rhosurf * custom_tstar);
204  q_star[lev]->setVal(custom_rhosurf * custom_qstar);
205  } else {
206  u_star[lev]->setVal(custom_ustar);
207  t_star[lev]->setVal(custom_tstar);
208  q_star[lev]->setVal(custom_qstar);
209  }
210  }
211 
212  if (m_update_k_rans) {
213  const bool use_ref_theta = (theta_ref > 0);
215  const Real l_inv_Cmu2 = inv_Cmu2;
216  const int klo = m_geom[lev].Domain().smallEnd(2);
217  IntVect ng = u_star[lev]->nGrowVect(); ng[2] = 0;
218 
219  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi)
220  {
221  Box gpbx = mfi.tilebox(IntVect(0),ng);
222 
223  if (gpbx.smallEnd(2) != klo) { continue; }
224 
225  gpbx.makeSlab(2,klo);
226 
227  auto cons_arr = cons_in.array(mfi);
228  const auto& u_star_arr = u_star[lev]->const_array(mfi);
229  const auto& t_star_arr = t_star[lev]->const_array(mfi);
230  const auto& dist_arr = walldist->const_array(mfi);
231 
232  ParallelFor(gpbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
233  {
234  Real rho = cons_arr(i,j,k,Rho_comp);
235  if (t_star_arr(i,j,0) < -1e-8) {
236  // Only destabilizing buoyancy flux affects the boundary k
237  // tstar < 0 ==> B > 0
238  Real B = -CONST_GRAV * l_inv_theta0 * u_star_arr(i,j,0) * t_star_arr(i,j,0);
239  if (!use_ref_theta) {
240  B *= cons_arr(i,j,k,Rho_comp) /
241  cons_arr(i,j,k,RhoTheta_comp);
242  }
243 
244  // Axell & Liungman 2001, Eqn. 16
245  cons_arr(i,j,k,RhoKE_comp) = rho * l_inv_Cmu2 *
246  std::pow(
247  u_star_arr(i,j,0) * u_star_arr(i,j,0) * u_star_arr(i,j,0)
248  + KAPPA * B * dist_arr(i,j,k),
249  two/three);
250  } else {
251  cons_arr(i,j,k,RhoKE_comp) = rho * l_inv_Cmu2 * u_star_arr(i,j,0) * u_star_arr(i,j,0);
252  }
253  });
254  }
255  }
256 
257  u_star[lev]->FillBoundary(m_geom[lev].periodicity());
258  t_star[lev]->FillBoundary(m_geom[lev].periodicity());
259  q_star[lev]->FillBoundary(m_geom[lev].periodicity());
260  olen[lev]->FillBoundary(m_geom[lev].periodicity());
261 }
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real two
Definition: ERF_Constants.H:8
const bool use_ref_theta
Definition: ERF_SetupDiff.H:16
const Real l_inv_theta0
Definition: ERF_SetupDiff.H:17
void compute_averages(const int &lev)
Definition: ERF_MOSTAverage.cpp:900
void update_surf_temp(const amrex::Real &time)
Definition: ERF_SurfaceLayer.H:595
void get_lsm_tsurf(const int &lev)
Definition: ERF_SurfaceLayer.cpp:1106
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:1063
void compute_fluxes(const int &lev, const int &max_iters, amrex::MultiFab &cons_in, const FluxIter &most_flux, bool is_land)
void fill_tsurf_with_sst_and_tsk(const int &lev, const amrex::Real &time)
Definition: ERF_SurfaceLayer.cpp:976
void compute_sfc_params_from_lsm_fluxes(const int &lev, amrex::MultiFab &cons_in)
Definition: ERF_SurfaceLayer.cpp:907
Definition: ERF_MOSTStress.H:189
Definition: ERF_MOSTStress.H:367
Definition: ERF_EBMOSTStress.H:12
Definition: ERF_MOSTStress.H:282
Definition: ERF_MOSTStress.H:437
Definition: ERF_MOSTStress.H:140
Definition: ERF_MOSTStress.H:622
Definition: ERF_MOSTStress.H:852
Definition: ERF_EBMOSTStress.H:194
Definition: ERF_MOSTStress.H:741
Definition: ERF_MOSTStress.H:960
Definition: ERF_MOSTStress.H:517
Definition: ERF_MOSTStress.H:1215
Definition: ERF_MOSTStress.H:1547
Definition: ERF_EBMOSTStress.H:61
Definition: ERF_MOSTStress.H:1385
Definition: ERF_MOSTStress.H:1665
Definition: ERF_MOSTStress.H:1081
Here is the call graph for this function:

◆ 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
615  {
616  m_ma.update_field_ptrs(lev, vars_old, Theta_prim, Qv_prim, Qr_prim);
617  }
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:279
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 
)
1148 {
1150  MYNNPBLH estimator;
1151  compute_pblh(lev, vars, z_phys_cc, estimator, moisture_indices);
1153  amrex::Error("YSU/MRF PBLH calc not implemented yet");
1154  }
1155 }
void compute_pblh(const int &lev, amrex::Vector< amrex::Vector< amrex::MultiFab >> &vars, amrex::MultiFab *z_phys_cc, const PBLHeightEstimator &est, const MoistureComponentIndices &moisture_indice)
Definition: ERF_PBLHeight.H:8

◆ update_sst_ptr()

void SurfaceLayer::update_sst_ptr ( const int  lev,
const int  itime,
amrex::MultiFab *  sst_ptr 
)
inline
670  {
671  m_sst_lev[lev][itime] = sst_ptr;
672  }

◆ update_surf_temp()

void SurfaceLayer::update_surf_temp ( const amrex::Real time)
inline
596  {
597  if (m_has_ocean_lsm_tsurf) {
598  return;
599  }
600  if (surf_heating_rate != 0) {
601  int nlevs = static_cast<int>(m_geom.size());
602  for (int lev = 0; lev < nlevs; lev++) {
603  t_surf[lev]->setVal(surf_temp + surf_heating_rate * time);
604  amrex::Print() << "Surface temp at t=" << time << ": "
605  << surf_temp + surf_heating_rate * time << std::endl;
606  }
607  }
608  }

◆ update_tsk_ptr()

void SurfaceLayer::update_tsk_ptr ( const int  lev,
const int  itime,
amrex::MultiFab *  tsk_ptr 
)
inline
674  {
675  m_tsk_lev[lev][itime] = tsk_ptr;
676  }

Member Data Documentation

◆ cnk_a

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

◆ cnk_visc

bool SurfaceLayer::cnk_visc {false}
private

◆ custom_qstar

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

◆ custom_rhosurf

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

◆ custom_tstar

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

◆ custom_ustar

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

◆ default_land_surf_moist

amrex::Real SurfaceLayer::default_land_surf_moist {zero}
private

◆ default_land_surf_temp

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

◆ depth

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

◆ flux_type

FluxCalcType SurfaceLayer::flux_type {FluxCalcType::MOENG}

◆ inv_Cmu2

amrex::Real SurfaceLayer::inv_Cmu2 = zero
private

◆ m_Cd

amrex::Real SurfaceLayer::m_Cd = zero
private

◆ m_Ch

amrex::Real SurfaceLayer::m_Ch = zero
private

◆ m_Cq

amrex::Real SurfaceLayer::m_Cq = zero
private

◆ m_eb_vec

amrex::Vector<const eb_*> SurfaceLayer::m_eb_vec
private

◆ m_eddyDiffs_lev

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

◆ m_final_low_time

amrex::Real SurfaceLayer::m_final_low_time
private

◆ m_geom

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

Referenced by update_surf_temp().

◆ m_has_lsm_fluxes

bool SurfaceLayer::m_has_lsm_fluxes = false
private

◆ m_has_lsm_tsurf

bool SurfaceLayer::m_has_lsm_tsurf = false
private

◆ m_has_ocean_lsm_tsurf

bool SurfaceLayer::m_has_ocean_lsm_tsurf = false
private

◆ m_Hwave_lev

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

◆ m_ignore_sst

bool SurfaceLayer::m_ignore_sst = false
private

◆ m_include_wstar

bool SurfaceLayer::m_include_wstar = false
private

◆ m_lmask_lev

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

◆ m_low_time_interval

amrex::Real SurfaceLayer::m_low_time_interval
private

◆ m_lsm_data_lev

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

◆ m_lsm_data_name

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

◆ m_lsm_flux_lev

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

◆ m_lsm_flux_name

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

◆ m_lsm_tsurf_indx

int SurfaceLayer::m_lsm_tsurf_indx = -1
private

◆ m_Lwave_lev

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

◆ m_ma

MOSTAverage SurfaceLayer::m_ma
private

◆ m_rotate

bool SurfaceLayer::m_rotate = false
private

◆ m_sst_lev

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

◆ m_start_low_time

amrex::Real SurfaceLayer::m_start_low_time
private

◆ m_terrain_type

TerrainType SurfaceLayer::m_terrain_type
private

◆ m_tsk_lev

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

◆ m_update_k_rans

bool SurfaceLayer::m_update_k_rans = false
private

◆ m_var_z0

bool SurfaceLayer::m_var_z0 {false}
private

◆ moist_type

MoistCalcType SurfaceLayer::moist_type {MoistCalcType::ADIABATIC}

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

◆ q_star

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

◆ q_surf

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

◆ rico_qsat_z0

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

◆ rico_theta_z0

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

◆ rough_type_land

RoughCalcType SurfaceLayer::rough_type_land {RoughCalcType::CONSTANT}

◆ 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 update_surf_temp().

◆ surf_moist

amrex::Real SurfaceLayer::surf_moist
private

◆ surf_moist_flux

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

◆ surf_temp

amrex::Real SurfaceLayer::surf_temp
private

Referenced by update_surf_temp().

◆ surf_temp_flux

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

◆ t_star

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

◆ t_surf

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

◆ theta_ref

amrex::Real SurfaceLayer::theta_ref = zero
private

◆ theta_type

ThetaCalcType SurfaceLayer::theta_type {ThetaCalcType::ADIABATIC}

◆ u_star

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

◆ use_moisture

bool SurfaceLayer::use_moisture
private

◆ w_star

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

◆ z0_const

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

◆ z_0

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

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