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

#include <ERF_NOAHMP.H>

Inheritance diagram for NOAHMP:
Collaboration diagram for NOAHMP:

Public Member Functions

 NOAHMP ()
 
virtual ~NOAHMP ()=default
 
void Define (SolverChoice &) override
 
void Init (const int &lev, const amrex::MultiFab &cons_in, const amrex::Geometry &geom, const amrex::Real &dt) override
 
void Advance_With_State (const int &lev, amrex::MultiFab &cons_in, amrex::MultiFab &xvel_in, amrex::MultiFab &yvel_in, amrex::MultiFab *hfx3_out, amrex::MultiFab *qfx3_out, const amrex::Real &elapsed_time, const amrex::Real &dt, const int &nstep) override
 
void Plot_Landfile (const int &nstep) override
 
amrex::MultiFab * Lsm_Data_Ptr (const int &varIdx) override
 
amrex::MultiFab * Lsm_Flux_Ptr (const int &varIdx) override
 
amrex::Geometry Lsm_Geom () override
 
int Lsm_Data_Size () override
 
int Lsm_Flux_Size () override
 
std::string Lsm_DataName (const int &varIdx) override
 
int Lsm_DataIndex (std::string varname) override
 
std::string Lsm_FluxName (const int &varIdx) override
 
int Lsm_FluxIndex (std::string varname) override
 
- Public Member Functions inherited from NullSurf
 NullSurf ()
 
virtual ~NullSurf ()=default
 
virtual void Advance (const amrex::Real &)
 
virtual void Update_Micro_Vars (amrex::MultiFab &)
 
virtual void Update_State_Vars (amrex::MultiFab &)
 
virtual void Copy_State_to_Micro (const amrex::MultiFab &)
 
virtual void Copy_Micro_to_State (amrex::MultiFab &)
 
virtual std::unordered_map< std::string, std::string > & Lsm_WRFInputNames ()
 

Private Types

using FabPtr = std::shared_ptr< amrex::MultiFab >
 

Private Attributes

int m_lsm_data_size = LsmData_NOAHMP::NumVars
 
int m_lsm_flux_size = LsmFlux_NOAHMP::NumVars
 
amrex::Vector< int > LsmDataMap
 
amrex::Vector< int > LsmFluxMap
 
amrex::Vector< std::string > LsmDataName
 
amrex::Vector< std::string > LsmFluxName
 
amrex::Geometry m_geom
 
amrex::Geometry m_lsm_geom
 
amrex::Real m_dt
 
int khi_lsm
 
int m_nz_lsm = 1
 
amrex::Real m_dz_lsm = one
 
amrex::Array< FabPtr, LsmData_NOAHMP::NumVarslsm_fab_data
 
amrex::Array< FabPtr, LsmFlux_NOAHMP::NumVarslsm_fab_flux
 
NoahmpIO_vector noahmpio_vect
 
amrex::Vector< std::unique_ptr< amrex::FArrayBox > > noahmp_input_tmp
 
amrex::Vector< std::unique_ptr< amrex::FArrayBox > > noahmp_output_tmp
 
int m_plot_int_1 = -1
 

Member Typedef Documentation

◆ FabPtr

using NOAHMP::FabPtr = std::shared_ptr<amrex::MultiFab>
private

Constructor & Destructor Documentation

◆ NOAHMP()

NOAHMP::NOAHMP ( )
inline
90 {}

◆ ~NOAHMP()

virtual NOAHMP::~NOAHMP ( )
virtualdefault

Member Function Documentation

◆ Advance_With_State()

void NOAHMP::Advance_With_State ( const int &  lev,
amrex::MultiFab &  cons_in,
amrex::MultiFab &  xvel_in,
amrex::MultiFab &  yvel_in,
amrex::MultiFab *  hfx3_out,
amrex::MultiFab *  qfx3_out,
const amrex::Real elapsed_time,
const amrex::Real dt,
const int &  nstep 
)
overridevirtual

Reimplemented from NullSurf.

241 {
242  // Verify we need to take another LSM step
243  Real NOAH_time = static_cast<Real>(noahmpio_vect[0].itimestep-1) * static_cast<Real>(noahmpio_vect[0].DTBL);
244  if (elapsed_time < NOAH_time) { return; }
245 
246  Box domain = m_geom.Domain();
247 
248  Print () << "Noah-MP driver started at time step: " << nstep+1 << std::endl;
249 
250  bool is_moist = (cons_in.nComp() > RhoQ1_comp);
251 
252  int klo = domain.smallEnd(2);
253 
254  // Loop over blocks to copy forcing data to Noahmp, drive the land model,
255  // and copy data back to ERF Multifabs.
256  int idb = 0;
257  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi, ++idb) {
258 
259  Box bx = mfi.tilebox();
260  Box gbx = mfi.tilebox(IntVect(0,0,0),IntVect(1,1,0));
261 
262  // Check if tile is at the lower boundary in lower z direction
263  if (bx.smallEnd(2) != klo) { continue; }
264 
265  bx.makeSlab(2,klo);
266  gbx.makeSlab(2,klo);
267 
268  // For limiting when populating ghost cells
269  int i_lo = bx.smallEnd(0); int i_hi = bx.bigEnd(0);
270  int j_lo = bx.smallEnd(1); int j_hi = bx.bigEnd(1);
271 
272  NoahmpIO_type* noahmpio = &noahmpio_vect[idb];
273 
274  const Array4<const Real>& U_PHY = xvel_in.const_array(mfi);
275  const Array4<const Real>& V_PHY = yvel_in.const_array(mfi);
276  const Array4<const Real>& CONS = cons_in.const_array(mfi);
277 
278  // Into NOAH-MP
279  const Array4<const Real>& SWDOWN = lsm_fab_data[LsmData_NOAHMP::sw_flux_dn]->const_array(mfi);
280  const Array4<const Real>& GLW = lsm_fab_data[LsmData_NOAHMP::lw_flux_dn]->const_array(mfi);
281  const Array4<const Real>& COSZEN = lsm_fab_data[LsmData_NOAHMP::cos_zenith_angle]->const_array(mfi);
282 
283  // Out of NOAH-MP
284  Array4<Real> TSK = lsm_fab_data[LsmData_NOAHMP::t_sfc]->array(mfi);
285  Array4<Real> EMISS = lsm_fab_data[LsmData_NOAHMP::sfc_emis]->array(mfi);
286  Array4<Real> ALBSFCDIR_VIS = lsm_fab_data[LsmData_NOAHMP::sfc_alb_dir_vis]->array(mfi);
287  Array4<Real> ALBSFCDIR_NIR = lsm_fab_data[LsmData_NOAHMP::sfc_alb_dir_nir]->array(mfi);
288  Array4<Real> ALBSFCDIF_VIS = lsm_fab_data[LsmData_NOAHMP::sfc_alb_dif_vis]->array(mfi);
289  Array4<Real> ALBSFCDIF_NIR = lsm_fab_data[LsmData_NOAHMP::sfc_alb_dif_nir]->array(mfi);
290 
291  // NOTE: Need to expose stresses and get stresses from NOAHMP
292  Array4<Real> q_flux_arr = lsm_fab_flux[LsmFlux_NOAHMP::q_flux]->array(mfi);
293  Array4<Real> t_flux_arr = lsm_fab_flux[LsmFlux_NOAHMP::t_flux]->array(mfi);
294  Array4<Real> tau13_arr = lsm_fab_flux[LsmFlux_NOAHMP::tau13]->array(mfi);
295  Array4<Real> tau23_arr = lsm_fab_flux[LsmFlux_NOAHMP::tau23]->array(mfi);
296 
297  // Use The_Pinned_Arena() for host-accessible memory that can be used with GPU
298  Array4<Real> noah_input_arr = noahmp_input_tmp[idb]->array();
299  Array4<Real> noah_output_arr = noahmp_output_tmp[idb]->array();
300 
301  // Copy forcing data from ERF to Noahmp.
302  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
303  {
304  Real qv = (is_moist) ? CONS(i,j,k,RhoQ1_comp)/CONS(i,j,k,Rho_comp) : zero;
305  noah_input_arr(i,j,0,NoahmpInputComp::u_phy) = myhalf*(U_PHY(i,j,k)+U_PHY(i+1,j,k));
306  noah_input_arr(i,j,0,NoahmpInputComp::v_phy) = myhalf*(V_PHY(i,j,k)+V_PHY(i ,j+1,k));
307  noah_input_arr(i,j,0,NoahmpInputComp::t_phy) = getTgivenRandRTh(CONS(i,j,k,Rho_comp),CONS(i,j,k,RhoTheta_comp),qv);
308  noah_input_arr(i,j,0,NoahmpInputComp::qv_curr) = qv;
309  noah_input_arr(i,j,0,NoahmpInputComp::p8w) = getPgivenRTh(CONS(i,j,k,RhoTheta_comp),qv);
310  noah_input_arr(i,j,0,NoahmpInputComp::swdown) = SWDOWN(i,j,0);
311  noah_input_arr(i,j,0,NoahmpInputComp::glw) = GLW(i,j,0);
312  noah_input_arr(i,j,0,NoahmpInputComp::coszen) = COSZEN(i,j,0);
313  });
314 
315  // Synchronize to ensure GPU kernel is complete before host access
316  Gpu::streamSynchronize();
317 
318  // Now on the host, copy data to NoahmpIO arrays
319  LoopOnCpu(bx, [&] (int i, int j, int ) noexcept
320  {
321  noahmpio->U_PHY(i,1,j) = noah_input_arr(i,j,0,NoahmpInputComp::u_phy);
322  noahmpio->V_PHY(i,1,j) = noah_input_arr(i,j,0,NoahmpInputComp::v_phy);
323  noahmpio->T_PHY(i,1,j) = noah_input_arr(i,j,0,NoahmpInputComp::t_phy);
324  noahmpio->QV_CURR(i,1,j) = noah_input_arr(i,j,0,NoahmpInputComp::qv_curr);
325  noahmpio->P8W(i,1,j) = noah_input_arr(i,j,0,NoahmpInputComp::p8w);
326  noahmpio->SWDOWN(i,j) = noah_input_arr(i,j,0,NoahmpInputComp::swdown);
327  noahmpio->GLW(i,j) = noah_input_arr(i,j,0,NoahmpInputComp::glw);
328  noahmpio->COSZEN(i,j) = noah_input_arr(i,j,0,NoahmpInputComp::coszen);
329  });
330 
331  // Call the noahmpio driver code. This runs the land model forcing for
332  // each object in noahmpio_vect that represent a block in the domain.
333  noahmpio->itimestep += 1;
334  noahmpio->DriverMain();
335 
336  // Copy results from NoahmpIO back to temporary arrays
337  LoopOnCpu(bx, [&] (int i, int j, int ) noexcept
338  {
339  noah_output_arr(i,j,0,NoahmpOutputComp::hfx) = noahmpio->HFX(i,j);
340  noah_output_arr(i,j,0,NoahmpOutputComp::lh) = noahmpio->LH(i,j);
341  noah_output_arr(i,j,0,NoahmpOutputComp::tau_ew) = noahmpio->TAU_EW(i,j);
342  noah_output_arr(i,j,0,NoahmpOutputComp::tau_ns) = noahmpio->TAU_NS(i,j);
343  noah_output_arr(i,j,0,NoahmpOutputComp::tsk) = noahmpio->TSK(i,j);
344  noah_output_arr(i,j,0,NoahmpOutputComp::emiss) = noahmpio->EMISS(i,j);
345  noah_output_arr(i,j,0,NoahmpOutputComp::albsfcdir_vis) = noahmpio->ALBSFCDIRXY(i,1,j);
346  noah_output_arr(i,j,0,NoahmpOutputComp::albsfcdir_nir) = noahmpio->ALBSFCDIRXY(i,2,j);
347  noah_output_arr(i,j,0,NoahmpOutputComp::albsfcdif_vis) = noahmpio->ALBSFCDIFXY(i,1,j);
348  noah_output_arr(i,j,0,NoahmpOutputComp::albsfcdif_nir) = noahmpio->ALBSFCDIFXY(i,2,j);
349  });
350 
351  // Copy forcing data from Noahmp to ERF
352  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
353  {
354  // Limit indices to the valid box. FillBoundary will pick these up below.
355  int ii = std::min(std::max(i,i_lo),i_hi);
356  int jj = std::min(std::max(j,j_lo),j_hi);
357 
358  // SurfaceLayer fluxes at CC.
359  // Noah-MP returns the -9999 fill value for cells it does NOT process
360  // (sea-ice / open-water points, which still have LANDMASK=1). Applying
361  // that as a flux gives -9999/(rho*Cp) ~ -7.6 K*m/s and crashes the
362  // lowest cell to ~200 K. Detect the fill and instead write the
363  // lsm_flux_undefined sentinel; the surface layer then falls back to
364  // the MOST flux for those cells (see ERF_SurfaceLayer.cpp).
365  // NOTE: tau13/tau23 are nodal in xz/yz; the 2D MFs have 1 ghost cell
366  // so the surface layer can average them.
367  Real hfx_lsm = noah_output_arr(ii,jj,0,NoahmpOutputComp::hfx);
368  if (hfx_lsm > Real(-9990.0)) {
369  t_flux_arr(i,j,k) = hfx_lsm/(CONS(ii,jj,k,Rho_comp)*Cp_d);
370  q_flux_arr(i,j,k) = noah_output_arr(ii,jj,0,NoahmpOutputComp::lh)/(CONS(ii,jj,k,Rho_comp)*L_v);
371  tau13_arr(i,j,k) = noah_output_arr(ii,jj,0,NoahmpOutputComp::tau_ew)/CONS(ii,jj,k,Rho_comp);
372  tau23_arr(i,j,k) = noah_output_arr(ii,jj,0,NoahmpOutputComp::tau_ns)/CONS(ii,jj,k,Rho_comp);
373  } else {
374  t_flux_arr(i,j,k) = lsm_flux_undefined;
375  q_flux_arr(i,j,k) = lsm_flux_undefined;
376  tau13_arr(i,j,k) = lsm_flux_undefined;
377  tau23_arr(i,j,k) = lsm_flux_undefined;
378  }
379 
380  // RRTMGP variables
381  TSK(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::tsk);
382  EMISS(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::emiss);
383  ALBSFCDIR_VIS(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::albsfcdir_vis);
384  ALBSFCDIR_NIR(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::albsfcdir_nir);
385  ALBSFCDIF_VIS(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::albsfcdif_vis);
386  ALBSFCDIF_NIR(i,j,0) = noah_output_arr(ii,jj,0,NoahmpOutputComp::albsfcdif_nir);
387  });
388  }
389 
390  // Fill the ghost cells
391  for (auto ivar = 0; ivar < LsmFlux_NOAHMP::NumVars; ++ivar) {
392  lsm_fab_flux[ivar]->FillBoundary(m_geom.periodicity());
393  }
394  Print () << "Noah-MP driver completed" << std::endl;
395 };
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:31
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
constexpr amrex::Real lsm_flux_undefined
Definition: ERF_Constants.H:22
constexpr amrex::Real L_v
Definition: ERF_Constants.H:35
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
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
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
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::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Array< FabPtr, LsmFlux_NOAHMP::NumVars > lsm_fab_flux
Definition: ERF_NOAHMP.H:241
amrex::Vector< std::unique_ptr< amrex::FArrayBox > > noahmp_output_tmp
Definition: ERF_NOAHMP.H:248
amrex::Vector< std::unique_ptr< amrex::FArrayBox > > noahmp_input_tmp
Definition: ERF_NOAHMP.H:247
NoahmpIO_vector noahmpio_vect
Definition: ERF_NOAHMP.H:244
amrex::Array< FabPtr, LsmData_NOAHMP::NumVars > lsm_fab_data
Definition: ERF_NOAHMP.H:238
amrex::Geometry m_geom
Definition: ERF_NOAHMP.H:220
@ sw_flux_dn
Definition: ERF_NOAHMP.H:32
@ sfc_emis
Definition: ERF_NOAHMP.H:26
@ sfc_alb_dir_vis
Definition: ERF_NOAHMP.H:27
@ sfc_alb_dif_nir
Definition: ERF_NOAHMP.H:30
@ lw_flux_dn
Definition: ERF_NOAHMP.H:37
@ cos_zenith_angle
Definition: ERF_NOAHMP.H:31
@ t_sfc
Definition: ERF_NOAHMP.H:25
@ sfc_alb_dir_nir
Definition: ERF_NOAHMP.H:28
@ sfc_alb_dif_vis
Definition: ERF_NOAHMP.H:29
@ t_flux
Definition: ERF_NOAHMP.H:45
@ tau13
Definition: ERF_NOAHMP.H:47
@ q_flux
Definition: ERF_NOAHMP.H:46
@ NumVars
Definition: ERF_NOAHMP.H:49
@ tau23
Definition: ERF_NOAHMP.H:48
@ qv
Definition: ERF_Kessler.H:29
@ v_phy
Definition: ERF_NOAHMP.H:57
@ glw
Definition: ERF_NOAHMP.H:62
@ qv_curr
Definition: ERF_NOAHMP.H:59
@ u_phy
Definition: ERF_NOAHMP.H:56
@ coszen
Definition: ERF_NOAHMP.H:63
@ t_phy
Definition: ERF_NOAHMP.H:58
@ swdown
Definition: ERF_NOAHMP.H:61
@ p8w
Definition: ERF_NOAHMP.H:60
@ hfx
Definition: ERF_NOAHMP.H:70
@ tau_ns
Definition: ERF_NOAHMP.H:73
@ tau_ew
Definition: ERF_NOAHMP.H:72
@ tsk
Definition: ERF_NOAHMP.H:74
@ albsfcdir_vis
Definition: ERF_NOAHMP.H:76
@ albsfcdif_nir
Definition: ERF_NOAHMP.H:79
@ emiss
Definition: ERF_NOAHMP.H:75
@ albsfcdif_vis
Definition: ERF_NOAHMP.H:78
@ albsfcdir_nir
Definition: ERF_NOAHMP.H:77
@ lh
Definition: ERF_NOAHMP.H:71
Here is the call graph for this function:

◆ Define()

void NOAHMP::Define ( SolverChoice )
inlineoverridevirtual

Reimplemented from NullSurf.

98  {
99  // NOTE: We should parse constants from sc here if needed,
100  }

◆ Init()

void NOAHMP::Init ( const int &  lev,
const amrex::MultiFab &  cons_in,
const amrex::Geometry &  geom,
const amrex::Real dt 
)
overridevirtual

Reimplemented from NullSurf.

21 {
22 
23  m_dt = dt;
24  m_geom = geom;
25 
26  Box domain = geom.Domain();
27  khi_lsm = domain.smallEnd(2) - 1;
28 
38  LsmDataName = {"t_sfc" , "sfc_emis" ,
39  "sfc_alb_dir_vis" , "sfc_alb_dir_nir" ,
40  "sfc_alb_dif_vis" , "sfc_alb_dif_nir" ,
41  "cos_zenith_angle" , "sw_flux_dn" ,
42  "sw_flux_dn_dir_vis", "sw_flux_dn_dir_nir",
43  "sw_flux_dn_dif_vis", "sw_flux_dn_dif_nir",
44  "lw_flux_dn" };
45 
46 
51  LsmFluxName = {"t_flux" , "q_flux" ,
52  "tau13" , "tau23" };
53 
54  ParmParse pp("erf");
55  pp.query("plot_int_1" , m_plot_int_1);
56 
57  // NOTE: All boxes in ba extend from zlo to zhi, so this transform is valid.
58  // If that were to change, the dm and new ba are no longer valid and
59  // direct copying between lsm data/flux vars cannot be done in a parfor.
60 
61  // Set 2D box array for lsm data
62  IntVect ng(1,1,0);
63  BoxArray ba = cons_in.boxArray();
64  DistributionMapping dm = cons_in.DistributionMap();
65  BoxList bl_lsm = ba.boxList();
66  for (auto& b : bl_lsm) { b.setRange(2,0); }
67  BoxArray ba_lsm(std::move(bl_lsm));
68 
69  // Set up lsm geometry
70  const RealBox& dom_rb = m_geom.ProbDomain();
71  const Real* dom_dx = m_geom.CellSize();
72  RealBox lsm_rb = dom_rb;
73  Real lsm_dx[AMREX_SPACEDIM] = {AMREX_D_DECL(dom_dx[0],dom_dx[1],m_dz_lsm)};
74  Real lsm_z_hi = dom_rb.lo(2);
75  Real lsm_z_lo = lsm_z_hi - Real(m_nz_lsm)*lsm_dx[2];
76  lsm_rb.setHi(2,lsm_z_hi); lsm_rb.setLo(2,lsm_z_lo);
77  m_lsm_geom.define( ba_lsm.minimalBox(), lsm_rb, m_geom.Coord(), m_geom.isPeriodic() );
78 
79  // Create the data
80  for (auto ivar = 0; ivar < LsmData_NOAHMP::NumVars; ++ivar) {
81  // State vars are CC
82  lsm_fab_data[ivar] = std::make_shared<MultiFab>(ba_lsm, dm, 1, ng);
83 
84  // NOTE: Radiation steps first so we set values
85  // to reasonable initialization for coupling
86  Real val_to_set = zero;
87  if (ivar == LsmData_NOAHMP::t_sfc) {
88  val_to_set = Real(300.0);
89  } else if (ivar == LsmData_NOAHMP::sfc_emis) {
90  val_to_set = Real(0.9);
91  } else if ( (ivar>=LsmData_NOAHMP::sfc_alb_dir_vis) &&
93  val_to_set = Real(0.06);
94  } else {
95  val_to_set = zero;
96  }
97  lsm_fab_data[ivar]->setVal(val_to_set);
98  }
99 
100  // Create the fluxes
101  for (auto ivar = 0; ivar < LsmFlux_NOAHMP::NumVars; ++ivar) {
102  // NOTE: Fluxes are CC with ghost cells for averaging
103  lsm_fab_flux[ivar] = std::make_shared<MultiFab>(ba_lsm, dm, 1, IntVect(1,1,0));
104  lsm_fab_flux[ivar]->setVal(0.);
105  }
106 
107  Print() << "Noah-MP initialization started" << std::endl;
108 
109  // Set noahmpio_vect to the size of local blocks (boxes)
110  noahmpio_vect.resize(cons_in.local_size(), lev);
111 
112  // Allocate pinned buffer space for all the boxes
113  noahmp_input_tmp.resize(cons_in.local_size());
114  noahmp_output_tmp.resize(cons_in.local_size());
115 
116  int klo = domain.smallEnd(2);
117 
118  // Iterate over multifab and noahmpio object together. Multifabs is
119  // used to extract size of blocks and set bounds for noahmpio objects.
120  int idb = 0;
121  for (MFIter mfi(cons_in); mfi.isValid(); ++mfi, ++idb) {
122 
123  // Get bounds for the tile
124  Box bx = mfi.tilebox();
125 
126  // Check if tile is at the lower boundary in lower z direction
127  if (bx.smallEnd(2) != klo) { continue; }
128 
129  // Make a slab
130  bx.makeSlab(2,klo);
131 
132  // Allocate pinned buffers for each box
133  noahmp_input_tmp[idb] = std::make_unique<FArrayBox>(bx, NoahmpInputComp::NumComps , The_Pinned_Arena());
134  noahmp_output_tmp[idb] = std::make_unique<FArrayBox>(bx, NoahmpOutputComp::NumComps, The_Pinned_Arena());
135 
136  // Get reference to the noahmpio object
137  NoahmpIO_type* noahmpio = &noahmpio_vect[idb];
138 
139  // Pass idb context to noahmpio
140  noahmpio->blkid = idb;
141 
142  // Pass level context to noahmpio
143  noahmpio->level = lev;
144 
145  // Initialize scalar values
146  noahmpio->ScalarInitDefault();
147 
148  // Store the rank of process for noahmp
149  noahmpio->rank = ParallelDescriptor::MyProc();
150 
151  // Store parallel communicator for noahmp
152  noahmpio->comm = MPI_Comm_c2f(ParallelDescriptor::Communicator());
153 
154  // Read namelist.erf file. This file contains
155  // noahmpio specific parameters and is read by
156  // the Fortran side of the implementation.
157  noahmpio->ReadNamelist();
158 
159  // Read the headers from the NetCDF land file. This is also
160  // implemented on the Fortran side of things currently.
161  noahmpio->ReadLandHeader();
162 
163  // Extract tile bounds and set them to their corresponding
164  // noahmpio variables. At present we will set all the variables
165  // corresponding to domain, memory, and tile to the same bounds.
166  // This will be changed later if we want to do special memory
167  // management for expensive use cases.
168  noahmpio->xstart = bx.smallEnd(0);
169  noahmpio->xend = bx.bigEnd(0);
170  noahmpio->ystart = bx.smallEnd(1);
171  noahmpio->yend = bx.bigEnd(1);
172 
173  // Domain bounds
174  noahmpio->ids = noahmpio->xstart;
175  noahmpio->ide = noahmpio->xend;
176  noahmpio->jds = noahmpio->ystart;
177  noahmpio->jde = noahmpio->yend;
178  noahmpio->kds = 1;
179  noahmpio->kde = 2;
180 
181  // Tile bounds
182  noahmpio->its = noahmpio->xstart;
183  noahmpio->ite = noahmpio->xend;
184  noahmpio->jts = noahmpio->ystart;
185  noahmpio->jte = noahmpio->yend;
186  noahmpio->kts = 1;
187  noahmpio->kte = 2;
188 
189  // Memory bounds
190  noahmpio->ims = noahmpio->xstart;
191  noahmpio->ime = noahmpio->xend;
192  noahmpio->jms = noahmpio->ystart;
193  noahmpio->jme = noahmpio->yend;
194  noahmpio->kms = 1;
195  noahmpio->kme = 2;
196 
197  // This procedure allocates memory in Fortran for IO variables
198  // using bounds that are set above and read from namelist.erf
199  // and headers from the NetCDF land file
200  noahmpio->VarInitDefault();
201 
202  // This reads NoahmpTable.TBL file which is another input file
203  // we need to set some IO variables.
204  noahmpio->ReadTable();
205 
206  // Read and initialize data from the NetCDF land file.
207  noahmpio->ReadLandMain();
208 
209  // Compute additional initial values that were not supplied
210  // by the NetCDF land file.
211  noahmpio->InitMain();
212 
213  // Write initial plotfile for land with the tag 0
214  Print() << "Noah-MP writing lnd.nc file at lev: " << lev << std::endl;
215  noahmpio->WriteLand(0);
216  }
218 
219  Print() << "Noah-MP initialization completed" << std::endl;
220 
221 };
ParmParse pp("prob")
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
int khi_lsm
Definition: ERF_NOAHMP.H:229
int m_plot_int_1
Definition: ERF_NOAHMP.H:251
amrex::Vector< std::string > LsmDataName
Definition: ERF_NOAHMP.H:214
amrex::Vector< int > LsmDataMap
Definition: ERF_NOAHMP.H:208
int m_lsm_data_size
Definition: ERF_NOAHMP.H:202
amrex::Geometry m_lsm_geom
Definition: ERF_NOAHMP.H:223
int m_nz_lsm
Definition: ERF_NOAHMP.H:232
int m_lsm_flux_size
Definition: ERF_NOAHMP.H:205
amrex::Vector< int > LsmFluxMap
Definition: ERF_NOAHMP.H:211
amrex::Vector< std::string > LsmFluxName
Definition: ERF_NOAHMP.H:217
amrex::Real m_dz_lsm
Definition: ERF_NOAHMP.H:235
amrex::Real m_dt
Definition: ERF_NOAHMP.H:226
@ sw_flux_dn_dif_nir
Definition: ERF_NOAHMP.H:36
@ sw_flux_dn_dir_vis
Definition: ERF_NOAHMP.H:33
@ sw_flux_dn_dir_nir
Definition: ERF_NOAHMP.H:34
@ NumVars
Definition: ERF_NOAHMP.H:38
@ sw_flux_dn_dif_vis
Definition: ERF_NOAHMP.H:35
@ ng
Definition: ERF_Morrison.H:48
@ NumComps
Definition: ERF_NOAHMP.H:64
@ NumComps
Definition: ERF_NOAHMP.H:80
Here is the call graph for this function:

◆ Lsm_Data_Ptr()

amrex::MultiFab* NOAHMP::Lsm_Data_Ptr ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullSurf.

127  {
128  int lsmIdx = LsmDataMap[varIdx];
129  AMREX_ALWAYS_ASSERT(lsmIdx < NOAHMP::m_lsm_data_size && lsmIdx>=0);
130  return lsm_fab_data[lsmIdx].get();
131  }
Here is the call graph for this function:

◆ Lsm_Data_Size()

int NOAHMP::Lsm_Data_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

148 { return NOAHMP::m_lsm_data_size; }

◆ Lsm_DataIndex()

int NOAHMP::Lsm_DataIndex ( std::string  varname)
inlineoverridevirtual

Reimplemented from NullSurf.

166  {
167  int varIdx = -1;
168  std::string lc_varname = amrex::toLower(varname);
169  for (int idx(0); idx<LsmData_NOAHMP::NumVars; ++idx) {
170  if (lc_varname == amrex::toLower(LsmDataName[idx])) {
171  varIdx = idx;
172  }
173  }
174  return varIdx;
175  }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int idx(int i, int j, int k, int nx, int ny)
Definition: ERF_InitForEnsemble.cpp:287
Here is the call graph for this function:

◆ Lsm_DataName()

std::string NOAHMP::Lsm_DataName ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullSurf.

157  {
158  int lsmIdx = LsmDataMap[varIdx];
159  AMREX_ALWAYS_ASSERT(lsmIdx < NOAHMP::m_lsm_data_size && lsmIdx>=0);
160  return LsmDataName[lsmIdx];
161  }
Here is the call graph for this function:

◆ Lsm_Flux_Ptr()

amrex::MultiFab* NOAHMP::Lsm_Flux_Ptr ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullSurf.

136  {
137  int lsmIdx = LsmFluxMap[varIdx];
138  AMREX_ALWAYS_ASSERT(lsmIdx < NOAHMP::m_lsm_flux_size && lsmIdx>=0);
139  return lsm_fab_flux[lsmIdx].get();
140  }
Here is the call graph for this function:

◆ Lsm_Flux_Size()

int NOAHMP::Lsm_Flux_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

152 { return NOAHMP::m_lsm_flux_size; }

◆ Lsm_FluxIndex()

int NOAHMP::Lsm_FluxIndex ( std::string  varname)
inlineoverridevirtual

Reimplemented from NullSurf.

189  {
190  int varIdx = -1;
191  std::string lc_varname = amrex::toLower(varname);
192  for (int idx(0); idx<LsmFlux_NOAHMP::NumVars; ++idx) {
193  if (lc_varname == amrex::toLower(LsmFluxName[idx])) {
194  varIdx = idx;
195  }
196  }
197  return varIdx;
198  }
Here is the call graph for this function:

◆ Lsm_FluxName()

std::string NOAHMP::Lsm_FluxName ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullSurf.

180  {
181  int lsmIdx = LsmFluxMap[varIdx];
182  AMREX_ALWAYS_ASSERT(lsmIdx < NOAHMP::m_lsm_flux_size && lsmIdx>=0);
183  return LsmFluxName[lsmIdx];
184  }
Here is the call graph for this function:

◆ Lsm_Geom()

amrex::Geometry NOAHMP::Lsm_Geom ( )
inlineoverridevirtual

Reimplemented from NullSurf.

144 { return m_lsm_geom; }

◆ Plot_Landfile()

void NOAHMP::Plot_Landfile ( const int &  nstep)
overridevirtual

Reimplemented from NullSurf.

225 {
226  for (NoahmpIO_type &noahmpio : noahmpio_vect) {
227  noahmpio.WriteLand(nstep);
228  }
229 }

Member Data Documentation

◆ khi_lsm

int NOAHMP::khi_lsm
private

◆ lsm_fab_data

amrex::Array<FabPtr, LsmData_NOAHMP::NumVars> NOAHMP::lsm_fab_data
private

Referenced by Lsm_Data_Ptr().

◆ lsm_fab_flux

amrex::Array<FabPtr, LsmFlux_NOAHMP::NumVars> NOAHMP::lsm_fab_flux
private

Referenced by Lsm_Flux_Ptr().

◆ LsmDataMap

amrex::Vector<int> NOAHMP::LsmDataMap
private

Referenced by Lsm_Data_Ptr(), and Lsm_DataName().

◆ LsmDataName

amrex::Vector<std::string> NOAHMP::LsmDataName
private

Referenced by Lsm_DataIndex(), and Lsm_DataName().

◆ LsmFluxMap

amrex::Vector<int> NOAHMP::LsmFluxMap
private

Referenced by Lsm_Flux_Ptr(), and Lsm_FluxName().

◆ LsmFluxName

amrex::Vector<std::string> NOAHMP::LsmFluxName
private

Referenced by Lsm_FluxIndex(), and Lsm_FluxName().

◆ m_dt

amrex::Real NOAHMP::m_dt
private

◆ m_dz_lsm

amrex::Real NOAHMP::m_dz_lsm = one
private

◆ m_geom

amrex::Geometry NOAHMP::m_geom
private

◆ m_lsm_data_size

int NOAHMP::m_lsm_data_size = LsmData_NOAHMP::NumVars
private

Referenced by Lsm_Data_Size().

◆ m_lsm_flux_size

int NOAHMP::m_lsm_flux_size = LsmFlux_NOAHMP::NumVars
private

Referenced by Lsm_Flux_Size().

◆ m_lsm_geom

amrex::Geometry NOAHMP::m_lsm_geom
private

Referenced by Lsm_Geom().

◆ m_nz_lsm

int NOAHMP::m_nz_lsm = 1
private

◆ m_plot_int_1

int NOAHMP::m_plot_int_1 = -1
private

◆ noahmp_input_tmp

amrex::Vector<std::unique_ptr<amrex::FArrayBox> > NOAHMP::noahmp_input_tmp
private

◆ noahmp_output_tmp

amrex::Vector<std::unique_ptr<amrex::FArrayBox> > NOAHMP::noahmp_output_tmp
private

◆ noahmpio_vect

NoahmpIO_vector NOAHMP::noahmpio_vect
private

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