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

#include <ERF_SLM.H>

Inheritance diagram for SLM:
Collaboration diagram for SLM:

Public Member Functions

 SLM ()
 
virtual ~SLM ()=default
 
void Define (SolverChoice &) override
 
void Init (const int &, 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 AdvanceSLM ()
 
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 Advance_With_State (const int &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab &, amrex::MultiFab *, amrex::MultiFab *, const amrex::Real &, const int &)
 
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 int Lsm_VarIndex (std::string)
 

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_SLM::NumVarslsm_fab_vars
 
amrex::Array< FabPtr, LsmVar_SLM::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 = 400.0
 
amrex::Real m_d_soil = m_k_soil / m_cp_soil
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ SLM()

SLM::SLM ( )
inline
32 {}

◆ ~SLM()

virtual SLM::~SLM ( )
virtualdefault

Member Function Documentation

◆ Advance()

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

Reimplemented from NullSurf.

56  {
57  m_dt = dt;
58  this->ComputeFluxes();
59  this->AdvanceSLM();
60  this->ComputeTsurf();
61  }
void ComputeTsurf()
Definition: ERF_SLM.cpp:64
amrex::Real m_dt
Definition: ERF_SLM.H:127
void ComputeFluxes()
Definition: ERF_SLM.cpp:83
void AdvanceSLM()
Definition: ERF_SLM.cpp:108
Here is the call graph for this function:

◆ AdvanceSLM()

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

Referenced by Advance().

Here is the caller graph for this function:

◆ ComputeFluxes()

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

Referenced by Advance().

Here is the caller graph for this function:

◆ ComputeTsurf()

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

Referenced by Advance().

Here is the caller graph for this function:

◆ Define()

void SLM::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 SLM::Init ( const int &  ,
const amrex::MultiFab &  cons_in,
const amrex::Geometry &  geom,
const amrex::Real &  dt 
)
overridevirtual

Reimplemented from NullSurf.

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

◆ Lsm_Data_Ptr()

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

Reimplemented from NullSurf.

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

◆ Lsm_Data_Size()

int SLM::Lsm_Data_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

99 { return SLM::m_lsm_size; }

◆ Lsm_Flux_Ptr()

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

Reimplemented from NullSurf.

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

◆ Lsm_Geom()

amrex::Geometry SLM::Lsm_Geom ( )
inlineoverridevirtual

Reimplemented from NullSurf.

95 { return m_lsm_geom; }

◆ Lsm_VarName()

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

Reimplemented from NullSurf.

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

Member Data Documentation

◆ khi_lsm

int SLM::khi_lsm
private

◆ lsm_fab_flux

amrex::Array<FabPtr, LsmVar_SLM::NumVars> SLM::lsm_fab_flux
private

Referenced by Lsm_Flux_Ptr().

◆ lsm_fab_vars

amrex::Array<FabPtr, LsmVar_SLM::NumVars> SLM::lsm_fab_vars
private

Referenced by Lsm_Data_Ptr().

◆ LsmVarMap

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

◆ LsmVarName

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

Referenced by Lsm_VarName().

◆ m_cp_soil

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

◆ m_d_soil

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

◆ m_dt

amrex::Real SLM::m_dt
private

Referenced by Advance().

◆ m_dz_lsm

amrex::Real SLM::m_dz_lsm = 0.1
private

◆ m_geom

amrex::Geometry SLM::m_geom
private

◆ m_k_soil

amrex::Real SLM::m_k_soil = 0.2
private

◆ m_lsm_geom

amrex::Geometry SLM::m_lsm_geom
private

Referenced by Lsm_Geom().

◆ m_lsm_size

int SLM::m_lsm_size = 1
private

Referenced by Lsm_Data_Size().

◆ m_nz_lsm

int SLM::m_nz_lsm = 30
private

◆ m_theta_dir

amrex::Real SLM::m_theta_dir = 400.0
private

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