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

#include <ERF_LandSurface.H>

Collaboration diagram for LandSurface:

Public Member Functions

 LandSurface ()
 
 ~LandSurface ()=default
 
void ReSize (const int &nlev)
 
template<class NewSurfModel >
void SetModel ()
 
void Define (const int &lev, SolverChoice &sc)
 
void Init (const int &lev, const amrex::MultiFab &cons_in, const amrex::Geometry &geom, const amrex::Real &dt_advance)
 
void Advance (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 &dt_advance, const int &nstep)
 
void Advance (const int &lev, const amrex::Real &dt_advance)
 
void Plot (const int &lev, const int &nstep)
 
void Update_Micro_Vars_Lev (const int &lev, amrex::MultiFab &cons_in)
 
void Update_State_Vars_Lev (const int &lev, amrex::MultiFab &cons_in)
 
amrex::MultiFab * Get_Data_Ptr (const int &lev, const int &varIdx)
 
amrex::MultiFab * Get_Flux_Ptr (const int &lev, const int &varIdx)
 
amrex::Geometry Get_Lsm_Geom (const int &lev)
 
int Get_Data_Size ()
 
int Get_Flux_Size ()
 
std::string Get_DataName (const int &varIdx)
 
int Get_DataIdx (const int &lev, std::string &varname)
 
std::string Get_FluxName (const int &varIdx)
 
int Get_FluxIdx (const int &lev, std::string &varname)
 
std::unordered_map< std::string, std::string > & Get_WRFInputNames ()
 
void Plot_Lsm_Data (amrex::Real time, const amrex::Vector< int > &level_steps, const amrex::Vector< amrex::IntVect > &ref_ratio)
 
template<class SurfModelType >
SurfModelType * get_model_lev (const int &lev)
 

Private Attributes

amrex::Vector< std::unique_ptr< NullSurf > > m_lsm_model
 
std::string plot_file_lsm {"plt_lsm_"}
 
amrex::Vector< amrex::Geometry > m_lsm_geom_lev
 
amrex::Vector< amrex::MultiFab > m_lsm_data_lev
 

Constructor & Destructor Documentation

◆ LandSurface()

LandSurface::LandSurface ( )
inline
19 { }

◆ ~LandSurface()

LandSurface::~LandSurface ( )
default

Member Function Documentation

◆ Advance() [1/2]

void LandSurface::Advance ( 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 dt_advance,
const int &  nstep 
)
inline
60  {
61  m_lsm_model[lev]->Advance_With_State(lev, cons_in, xvel_in, yvel_in, hfx3_out, qfx3_out, dt_advance, nstep);
62  }
amrex::Vector< std::unique_ptr< NullSurf > > m_lsm_model
Definition: ERF_LandSurface.H:161

◆ Advance() [2/2]

void LandSurface::Advance ( const int &  lev,
const amrex::Real dt_advance 
)
inline
66  {
67  m_lsm_model[lev]->Advance(dt_advance);
68  }

◆ Define()

void LandSurface::Define ( const int &  lev,
SolverChoice sc 
)
inline
38  {
39  m_lsm_model[lev]->Define(sc);
40  }

◆ Get_Data_Ptr()

amrex::MultiFab* LandSurface::Get_Data_Ptr ( const int &  lev,
const int &  varIdx 
)
inline
89 { return m_lsm_model[lev]->Lsm_Data_Ptr(varIdx); }

Referenced by Plot_Lsm_Data().

Here is the caller graph for this function:

◆ Get_Data_Size()

int LandSurface::Get_Data_Size ( )
inline
98 { return m_lsm_model[0]->Lsm_Data_Size(); }

Referenced by Plot_Lsm_Data().

Here is the caller graph for this function:

◆ Get_DataIdx()

int LandSurface::Get_DataIdx ( const int &  lev,
std::string &  varname 
)
inline
107 { return m_lsm_model[lev]->Lsm_DataIndex(varname); }

◆ Get_DataName()

std::string LandSurface::Get_DataName ( const int &  varIdx)
inline
104 { return m_lsm_model[0]->Lsm_DataName(varIdx); }

Referenced by Plot_Lsm_Data().

Here is the caller graph for this function:

◆ Get_Flux_Ptr()

amrex::MultiFab* LandSurface::Get_Flux_Ptr ( const int &  lev,
const int &  varIdx 
)
inline
92 { return m_lsm_model[lev]->Lsm_Flux_Ptr(varIdx); }

◆ Get_Flux_Size()

int LandSurface::Get_Flux_Size ( )
inline
101 { return m_lsm_model[0]->Lsm_Flux_Size(); }

◆ Get_FluxIdx()

int LandSurface::Get_FluxIdx ( const int &  lev,
std::string &  varname 
)
inline
113 { return m_lsm_model[lev]->Lsm_FluxIndex(varname); }

◆ Get_FluxName()

std::string LandSurface::Get_FluxName ( const int &  varIdx)
inline
110 { return m_lsm_model[0]->Lsm_FluxName(varIdx); }

◆ Get_Lsm_Geom()

amrex::Geometry LandSurface::Get_Lsm_Geom ( const int &  lev)
inline
95 { return m_lsm_model[lev]->Lsm_Geom(); }

Referenced by Plot_Lsm_Data().

Here is the caller graph for this function:

◆ get_model_lev()

template<class SurfModelType >
SurfModelType* LandSurface::get_model_lev ( const int &  lev)
inline
155  {
156  return static_cast<SurfModelType*>(m_lsm_model[lev].get());
157  }

◆ Get_WRFInputNames()

std::unordered_map<std::string,std::string>& LandSurface::Get_WRFInputNames ( )
inline
116 { return m_lsm_model[0]->Lsm_WRFInputNames(); }

◆ Init()

void LandSurface::Init ( const int &  lev,
const amrex::MultiFab &  cons_in,
const amrex::Geometry &  geom,
const amrex::Real dt_advance 
)
inline
47  {
48  m_lsm_model[lev]->Init(lev, cons_in, geom, dt_advance);
49  }

◆ Plot()

void LandSurface::Plot ( const int &  lev,
const int &  nstep 
)
inline
72  {
73  m_lsm_model[lev]->Plot_Landfile(nstep);
74  }

◆ Plot_Lsm_Data()

void LandSurface::Plot_Lsm_Data ( amrex::Real  time,
const amrex::Vector< int > &  level_steps,
const amrex::Vector< amrex::IntVect > &  ref_ratio 
)
inline
122  {
123  int nlev = m_lsm_model.size();
124  int nvar = this->Get_Data_Size();
125  amrex::Vector<std::string> varnames;
126 
127  // Only write if we have valid pointers
128  if (this->Get_Data_Ptr(0,0)) {
129  varnames.resize(nvar);
130  m_lsm_geom_lev.resize(nlev);
131  m_lsm_data_lev.resize(nlev);
132  std::string plotfilename = amrex::Concatenate(plot_file_lsm, level_steps[0], 5);
133  for (int lev(0); lev<nlev; ++lev) {
134  m_lsm_geom_lev[lev] = this->Get_Lsm_Geom(lev);
135 
136  amrex::MultiFab* mf_lsm = this->Get_Data_Ptr(lev,0);
137  amrex::IntVect ng(0,0,0);
138  amrex::BoxArray ba = mf_lsm->boxArray();
139  amrex::DistributionMapping dm = mf_lsm->DistributionMap();
140  m_lsm_data_lev[lev].define(ba, dm, nvar, ng);
141  for (int n(0); n<nvar; ++n) {
142  mf_lsm = this->Get_Data_Ptr(lev,n);
143  amrex::MultiFab::Copy(m_lsm_data_lev[lev],*(mf_lsm),0,n,1,0);
144  if (lev==0) varnames[n] = this->Get_DataName(n);
145  }
146  }
147  WriteMultiLevelPlotfile (plotfilename, nlev, GetVecOfConstPtrs(m_lsm_data_lev),
148  varnames, m_lsm_geom_lev, time, level_steps, ref_ratio);
149  m_lsm_geom_lev.clear();
150  m_lsm_data_lev.clear();
151  }
152  }
amrex::Vector< amrex::Geometry > m_lsm_geom_lev
Definition: ERF_LandSurface.H:167
int Get_Data_Size()
Definition: ERF_LandSurface.H:98
amrex::Geometry Get_Lsm_Geom(const int &lev)
Definition: ERF_LandSurface.H:95
amrex::Vector< amrex::MultiFab > m_lsm_data_lev
Definition: ERF_LandSurface.H:170
std::string Get_DataName(const int &varIdx)
Definition: ERF_LandSurface.H:104
amrex::MultiFab * Get_Data_Ptr(const int &lev, const int &varIdx)
Definition: ERF_LandSurface.H:89
std::string plot_file_lsm
Definition: ERF_LandSurface.H:164
@ ng
Definition: ERF_Morrison.H:48
Here is the call graph for this function:

◆ ReSize()

void LandSurface::ReSize ( const int &  nlev)
inline
24 { m_lsm_model.resize(nlev); }

◆ SetModel()

template<class NewSurfModel >
void LandSurface::SetModel ( )
inline
29  {
30  for (int lev(0); lev<m_lsm_model.size(); ++lev) {
31  m_lsm_model[lev] = std::make_unique<NewSurfModel>();
32  }
33  }

◆ Update_Micro_Vars_Lev()

void LandSurface::Update_Micro_Vars_Lev ( const int &  lev,
amrex::MultiFab &  cons_in 
)
inline
78  {
79  m_lsm_model[lev]->Update_Micro_Vars(cons_in);
80  }

◆ Update_State_Vars_Lev()

void LandSurface::Update_State_Vars_Lev ( const int &  lev,
amrex::MultiFab &  cons_in 
)
inline
84  {
85  m_lsm_model[lev]->Update_State_Vars(cons_in);
86  }

Member Data Documentation

◆ m_lsm_data_lev

amrex::Vector<amrex::MultiFab> LandSurface::m_lsm_data_lev
private

Referenced by Plot_Lsm_Data().

◆ m_lsm_geom_lev

amrex::Vector<amrex::Geometry> LandSurface::m_lsm_geom_lev
private

Referenced by Plot_Lsm_Data().

◆ m_lsm_model

◆ plot_file_lsm

std::string LandSurface::plot_file_lsm {"plt_lsm_"}
private

Referenced by Plot_Lsm_Data().


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