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

#include <ERF_MM5.H>

Inheritance diagram for MM5:
Collaboration diagram for MM5:

Public Member Functions

 MM5 ()
 
virtual ~MM5 ()=default
 
void Define (SolverChoice &) override
 
void Init (const amrex::MultiFab &cons_in, const amrex::Geometry &geom, const amrex::Real &dt) override
 
void Advance (const amrex::Real &dt) override
 
void ComputeTsurf ()
 
void ComputeFluxes ()
 
void AdvanceMM5 ()
 
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
 
std::string Lsm_VarName (const int &varIdx) override
 
- Public Member Functions inherited from NullSurf
 NullSurf ()
 
virtual ~NullSurf ()=default
 
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 &)
 

Private Types

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

Private Attributes

int m_lsm_size = 1
 
amrex::Vector< int > LsmVarMap
 
amrex::Vector< std::string > LsmVarName
 
amrex::Geometry m_geom
 
amrex::Geometry m_lsm_geom
 
amrex::Real m_dt
 
int khi_lsm
 
amrex::Array< FabPtr, LsmVar_MM5::NumVarslsm_fab_vars
 
amrex::Array< FabPtr, LsmVar_MM5::NumVarslsm_fab_flux
 
int m_nz_lsm = 30
 
amrex::Real m_dz_lsm = 0.1
 
amrex::Real m_cp_soil = 1.26e6
 
amrex::Real m_k_soil = 0.2
 
amrex::Real m_theta_dir = 300.0
 
amrex::Real m_d_soil = m_k_soil / m_cp_soil
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ MM5()

MM5::MM5 ( )
inline
32 {}

◆ ~MM5()

virtual MM5::~MM5 ( )
virtualdefault

Member Function Documentation

◆ Advance()

void MM5::Advance ( const amrex::Real &  dt)
inlineoverridevirtual

Reimplemented from NullSurf.

55  {
56  m_dt = dt;
57  this->ComputeFluxes();
58  this->AdvanceMM5();
59  this->ComputeTsurf();
60  }
void ComputeFluxes()
Definition: ERF_MM5.cpp:82
amrex::Real m_dt
Definition: ERF_MM5.H:126
void AdvanceMM5()
Definition: ERF_MM5.cpp:107
void ComputeTsurf()
Definition: ERF_MM5.cpp:63
Here is the call graph for this function:

◆ AdvanceMM5()

void MM5::AdvanceMM5 ( )
108 {
109  // Expose for GPU copy
110  Real dt = m_dt;
111  Real dzInv = m_lsm_geom.InvCellSize(2);
112 
113  for ( MFIter mfi(*(lsm_fab_vars[LsmVar_MM5::theta])); mfi.isValid(); ++mfi) {
114  auto box3d = mfi.tilebox();
115 
116  auto theta_array = lsm_fab_vars[LsmVar_MM5::theta]->array(mfi);
117  auto theta_flux = lsm_fab_flux[LsmVar_MM5::theta]->array(mfi);
118 
119  ParallelFor( box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
120  {
121  theta_array(i,j,k) += dt * ( theta_flux(i,j,k+1) - theta_flux(i,j,k) ) * dzInv;
122  });
123  }
124 }
amrex::Array< FabPtr, LsmVar_MM5::NumVars > lsm_fab_vars
Definition: ERF_MM5.H:138
amrex::Geometry m_lsm_geom
Definition: ERF_MM5.H:123
amrex::Array< FabPtr, LsmVar_MM5::NumVars > lsm_fab_flux
Definition: ERF_MM5.H:141
@ theta
Definition: ERF_MM5.H:20

Referenced by Advance().

Here is the caller graph for this function:

◆ ComputeFluxes()

void MM5::ComputeFluxes ( )
83 {
84  // Expose for GPU copy
85  int khi = khi_lsm;
86  Real Dsoil = m_d_soil;
87  Real dzInv = m_lsm_geom.InvCellSize(2);
88 
89  for ( MFIter mfi(*(lsm_fab_flux[LsmVar_MM5::theta])); mfi.isValid(); ++mfi) {
90  auto box3d = mfi.tilebox();
91 
92  // Do not overwrite the flux at the top (comes from MOST BC)
93  if (box3d.bigEnd(2) == khi+1) box3d.setBig(2,khi);
94 
95  auto theta_array = lsm_fab_vars[LsmVar_MM5::theta]->array(mfi);
96  auto theta_flux = lsm_fab_flux[LsmVar_MM5::theta]->array(mfi);
97 
98  ParallelFor( box3d, [=] AMREX_GPU_DEVICE (int i, int j, int k)
99  {
100  theta_flux(i,j,k) = Dsoil * ( theta_array(i,j,k) - theta_array(i,j,k-1) ) * dzInv;
101  });
102  }
103 }
amrex::Real m_d_soil
Definition: ERF_MM5.H:161
int khi_lsm
Definition: ERF_MM5.H:129

Referenced by Advance().

Here is the caller graph for this function:

◆ ComputeTsurf()

void MM5::ComputeTsurf ( )
64 {
65  // Expose for GPU copy
66  int khi = khi_lsm;
67 
68  for ( MFIter mfi(*(lsm_fab_vars[LsmVar_MM5::theta])); mfi.isValid(); ++mfi) {
69  auto box2d = mfi.tilebox(); box2d.makeSlab(2,khi);
70 
71  auto theta_array = lsm_fab_vars[LsmVar_MM5::theta]->array(mfi);
72 
73  ParallelFor( box2d, [=] AMREX_GPU_DEVICE (int i, int j, int )
74  {
75  theta_array(i,j,khi+1) = 1.5*theta_array(i,j,khi) - 0.5*theta_array(i,j,khi-1);
76  });
77  }
78 }

Referenced by Advance().

Here is the caller graph for this function:

◆ Define()

void MM5::Define ( SolverChoice )
inlineoverridevirtual

Reimplemented from NullSurf.

40  {
41  // NOTE: We should parse things from sc here,
42  // but they are hard coded because this
43  // is a demonstration for now.
44  }

◆ Init()

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

Reimplemented from NullSurf.

10 {
11  m_dt = dt;
12  m_geom = geom;
13 
14  Box domain = geom.Domain();
15  khi_lsm = domain.smallEnd(2) - 1;
16 
17  LsmVarMap.resize(m_lsm_size);
19 
20  LsmVarName.resize(m_lsm_size);
21  LsmVarName = {"theta"};
22 
23  // NOTE: All boxes in ba extend from zlo to zhi, so this transform is valid.
24  // If that were to change, the dm and new ba are no longer valid and
25  // direct copying between lsm data/flux vars cannot be done in a parfor.
26 
27  // Set box array for lsm data
28  IntVect ng(0,0,1);
29  BoxArray ba = cons_in.boxArray();
30  DistributionMapping dm = cons_in.DistributionMap();
31  BoxList bl_lsm = ba.boxList();
32  for (auto& b : bl_lsm) {
33  b.setBig(2,khi_lsm); // First point below the surface
34  b.setSmall(2,khi_lsm - m_nz_lsm + 1); // Last point below the surface
35  }
36  BoxArray ba_lsm(std::move(bl_lsm));
37 
38  // Set up lsm geometry
39  const RealBox& dom_rb = m_geom.ProbDomain();
40  const Real* dom_dx = m_geom.CellSize();
41  RealBox lsm_rb = dom_rb;
42  Real lsm_dx[AMREX_SPACEDIM] = {AMREX_D_DECL(dom_dx[0],dom_dx[1],m_dz_lsm)};
43  Real lsm_z_hi = dom_rb.lo(2);
44  Real lsm_z_lo = lsm_z_hi - Real(m_nz_lsm)*lsm_dx[2];
45  lsm_rb.setHi(2,lsm_z_hi); lsm_rb.setLo(2,lsm_z_lo);
46  m_lsm_geom.define( ba_lsm.minimalBox(), lsm_rb, m_geom.Coord(), m_geom.isPeriodic() );
47 
48  // Create the data and fluxes
49  for (auto ivar = 0; ivar < LsmVar_MM5::NumVars; ++ivar) {
50  // State vars are CC
51  Real theta_0 = m_theta_dir;
52  lsm_fab_vars[ivar] = std::make_shared<MultiFab>(ba_lsm, dm, 1, ng);
53  lsm_fab_vars[ivar]->setVal(theta_0);
54 
55  // Fluxes are nodal in z
56  lsm_fab_flux[ivar] = std::make_shared<MultiFab>(convert(ba_lsm, IntVect(0,0,1)), dm, 1, IntVect(0,0,0));
57  lsm_fab_flux[ivar]->setVal(0.);
58  }
59 }
int m_lsm_size
Definition: ERF_MM5.H:111
amrex::Vector< int > LsmVarMap
Definition: ERF_MM5.H:114
amrex::Real m_dz_lsm
Definition: ERF_MM5.H:149
amrex::Geometry m_geom
Definition: ERF_MM5.H:120
int m_nz_lsm
Definition: ERF_MM5.H:146
amrex::Real m_theta_dir
Definition: ERF_MM5.H:158
amrex::Vector< std::string > LsmVarName
Definition: ERF_MM5.H:117
@ NumVars
Definition: ERF_MM5.H:21

◆ Lsm_Data_Ptr()

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

Reimplemented from NullSurf.

77  {
78  int lsmIdx = LsmVarMap[varIdx];
79  AMREX_ALWAYS_ASSERT(lsmIdx < MM5::m_lsm_size && lsmIdx>=0);
80  return lsm_fab_vars[lsmIdx].get();
81  }

◆ Lsm_Data_Size()

int MM5::Lsm_Data_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

98 { return MM5::m_lsm_size; }

◆ Lsm_Flux_Ptr()

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

Reimplemented from NullSurf.

86  {
87  int lsmIdx = LsmVarMap[varIdx];
88  AMREX_ALWAYS_ASSERT(lsmIdx < MM5::m_lsm_size && lsmIdx>=0);
89  return lsm_fab_flux[lsmIdx].get();
90  }

◆ Lsm_Geom()

amrex::Geometry MM5::Lsm_Geom ( )
inlineoverridevirtual

Reimplemented from NullSurf.

94 { return m_lsm_geom; }

◆ Lsm_VarName()

std::string MM5::Lsm_VarName ( const int &  varIdx)
inlineoverridevirtual

Reimplemented from NullSurf.

103  {
104  int lsmIdx = LsmVarMap[varIdx];
105  AMREX_ALWAYS_ASSERT(lsmIdx < MM5::m_lsm_size && lsmIdx>=0);
106  return LsmVarName[lsmIdx];
107  }

Member Data Documentation

◆ khi_lsm

int MM5::khi_lsm
private

◆ lsm_fab_flux

amrex::Array<FabPtr, LsmVar_MM5::NumVars> MM5::lsm_fab_flux
private

Referenced by Lsm_Flux_Ptr().

◆ lsm_fab_vars

amrex::Array<FabPtr, LsmVar_MM5::NumVars> MM5::lsm_fab_vars
private

Referenced by Lsm_Data_Ptr().

◆ LsmVarMap

amrex::Vector<int> MM5::LsmVarMap
private

◆ LsmVarName

amrex::Vector<std::string> MM5::LsmVarName
private

Referenced by Lsm_VarName().

◆ m_cp_soil

amrex::Real MM5::m_cp_soil = 1.26e6
private

◆ m_d_soil

amrex::Real MM5::m_d_soil = m_k_soil / m_cp_soil
private

◆ m_dt

amrex::Real MM5::m_dt
private

Referenced by Advance().

◆ m_dz_lsm

amrex::Real MM5::m_dz_lsm = 0.1
private

◆ m_geom

amrex::Geometry MM5::m_geom
private

◆ m_k_soil

amrex::Real MM5::m_k_soil = 0.2
private

◆ m_lsm_geom

amrex::Geometry MM5::m_lsm_geom
private

Referenced by Lsm_Geom().

◆ m_lsm_size

int MM5::m_lsm_size = 1
private

Referenced by Lsm_Data_Size().

◆ m_nz_lsm

int MM5::m_nz_lsm = 30
private

◆ m_theta_dir

amrex::Real MM5::m_theta_dir = 300.0
private

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