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

#include <ERF_OceanSurf.H>

Inheritance diagram for OceanSurf:
Collaboration diagram for OceanSurf:

Public Member Functions

 OceanSurf ()=default
 
 ~OceanSurf () override=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
 
amrex::MultiFab * Lsm_Data_Ptr (const int &varIdx) override
 
amrex::MultiFab * Lsm_Flux_Ptr (const int &) override
 
amrex::Geometry Lsm_Geom () override
 
int Lsm_Data_Size () override
 
int Lsm_Flux_Size () override
 
std::string Lsm_DataName (const int &varIdx) override
 
std::string Lsm_FluxName (const int &) override
 
int Lsm_DataIndex (std::string varname) override
 
int Lsm_FluxIndex (std::string) override
 
void SetOceanState (const amrex::Vector< amrex::MultiFab * > *ocean_state)
 
const amrex::Vector< amrex::MultiFab * > * OceanState () const
 
- 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 amrex::Real &, const int &)
 
virtual void Plot_Landfile (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 std::unordered_map< std::string, std::string > & Lsm_WRFInputNames ()
 

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::Vector< FabPtrlsm_fab_data
 
amrex::Geometry m_lsm_geom
 
amrex::Real m_dt {0.0}
 
int khi_lsm {0}
 
amrex::Real m_default_tsurf {amrex::Real(300.0)}
 
const amrex::Vector< amrex::MultiFab * > * m_ocean_state = nullptr
 

Member Typedef Documentation

◆ FabPtr

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

Constructor & Destructor Documentation

◆ OceanSurf()

OceanSurf::OceanSurf ( )
default

◆ ~OceanSurf()

OceanSurf::~OceanSurf ( )
overridedefault

Member Function Documentation

◆ Advance()

void OceanSurf::Advance ( const amrex::Real dt)
inlineoverridevirtual

Reimplemented from NullSurf.

72  {
73  // OceanSurf currently does not evolve internal state here; keep the
74  // latest timestep for bookkeeping so this stays aligned with other LSMs.
75  m_dt = dt;
76  }
amrex::Real m_dt
Definition: ERF_OceanSurf.H:152

◆ Define()

void OceanSurf::Define ( SolverChoice )
inlineoverridevirtual

Reimplemented from NullSurf.

26  {
27  amrex::ParmParse pp("erf");
29  pp.query("most.surf_temp", m_default_tsurf);
30  }
ParmParse pp("prob")
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Real m_default_tsurf
Definition: ERF_OceanSurf.H:154
Here is the call graph for this function:

◆ Init()

void OceanSurf::Init ( const int &  ,
const amrex::MultiFab &  cons_in,
const amrex::Geometry &  geom,
const amrex::Real dt 
)
inlineoverridevirtual

Reimplemented from NullSurf.

37  {
38  m_dt = dt;
39  m_lsm_geom = geom;
40 
41  amrex::Box domain = geom.Domain();
42  khi_lsm = domain.smallEnd(2);
43 
44  LsmVarMap.resize(m_lsm_size);
46 
47  LsmVarName.resize(m_lsm_size);
48  LsmVarName = {"t_surf"};
49 
50  amrex::IntVect ng(1,1,0);
51  amrex::BoxArray ba = cons_in.boxArray();
52  amrex::DistributionMapping dm = cons_in.DistributionMap();
53  amrex::BoxList bl_lsm = ba.boxList();
54  for (auto& b : bl_lsm) {
55  b.setRange(2,0);
56  }
57  amrex::BoxArray ba_lsm(std::move(bl_lsm));
58 
59  const amrex::RealBox& dom_rb = geom.ProbDomain();
60  amrex::RealBox lsm_rb = dom_rb;
61  lsm_rb.setHi(2, dom_rb.hi(2));
62  lsm_rb.setLo(2, dom_rb.lo(2));
63  m_lsm_geom.define(ba_lsm.minimalBox(), lsm_rb, geom.Coord(), geom.isPeriodic());
64 
65  lsm_fab_data.resize(m_lsm_size);
66  lsm_fab_data[0] = std::make_shared<amrex::MultiFab>(ba_lsm, dm, 1, ng);
67  lsm_fab_data[0]->setVal(m_default_tsurf);
68  }
amrex::Vector< FabPtr > lsm_fab_data
Definition: ERF_OceanSurf.H:150
amrex::Vector< std::string > LsmVarName
Definition: ERF_OceanSurf.H:149
int m_lsm_size
Definition: ERF_OceanSurf.H:147
amrex::Geometry m_lsm_geom
Definition: ERF_OceanSurf.H:151
int khi_lsm
Definition: ERF_OceanSurf.H:153
amrex::Vector< int > LsmVarMap
Definition: ERF_OceanSurf.H:148
@ t_surf
Definition: ERF_OceanSurf.H:14
@ ng
Definition: ERF_Morrison.H:48

◆ Lsm_Data_Ptr()

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

Reimplemented from NullSurf.

80  {
81  int lsmIdx = LsmVarMap[varIdx];
82  AMREX_ALWAYS_ASSERT(lsmIdx < OceanSurf::m_lsm_size && lsmIdx>=0);
83  return lsm_fab_data[lsmIdx].get();
84  }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
Here is the call graph for this function:

◆ Lsm_Data_Size()

int OceanSurf::Lsm_Data_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

96 { return OceanSurf::m_lsm_size; }

◆ Lsm_DataIndex()

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

Reimplemented from NullSurf.

117  {
118  std::string lc_varname = amrex::toLower(varname);
119  if (lc_varname == "t_surf" || lc_varname == "theta") {
120  return 0;
121  }
122  return -1;
123  }

◆ Lsm_DataName()

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

Reimplemented from NullSurf.

103  {
104  int lsmIdx = LsmVarMap[varIdx];
105  AMREX_ALWAYS_ASSERT(lsmIdx < OceanSurf::m_lsm_size && lsmIdx>=0);
106  return LsmVarName[lsmIdx];
107  }
Here is the call graph for this function:

◆ Lsm_Flux_Ptr()

amrex::MultiFab* OceanSurf::Lsm_Flux_Ptr ( const int &  )
inlineoverridevirtual

Reimplemented from NullSurf.

88  {
89  return nullptr;
90  }

◆ Lsm_Flux_Size()

int OceanSurf::Lsm_Flux_Size ( )
inlineoverridevirtual

Reimplemented from NullSurf.

99 { return 0; }

◆ Lsm_FluxIndex()

int OceanSurf::Lsm_FluxIndex ( std::string  )
inlineoverridevirtual

Reimplemented from NullSurf.

127  {
128  return -1;
129  }

◆ Lsm_FluxName()

std::string OceanSurf::Lsm_FluxName ( const int &  )
inlineoverridevirtual

Reimplemented from NullSurf.

111  {
112  return "";
113  }

◆ Lsm_Geom()

amrex::Geometry OceanSurf::Lsm_Geom ( )
inlineoverridevirtual

Reimplemented from NullSurf.

93 { return m_lsm_geom; }

◆ OceanState()

const amrex::Vector<amrex::MultiFab*>* OceanSurf::OceanState ( ) const
inline
140  {
141  return m_ocean_state;
142  }
const amrex::Vector< amrex::MultiFab * > * m_ocean_state
Definition: ERF_OceanSurf.H:155

◆ SetOceanState()

void OceanSurf::SetOceanState ( const amrex::Vector< amrex::MultiFab * > *  ocean_state)
inline
133  {
134  // This pointer is currently packed and filled by ERF_Coupling.cpp.
135  m_ocean_state = ocean_state;
136  }

Member Data Documentation

◆ khi_lsm

int OceanSurf::khi_lsm {0}
private

Referenced by Init().

◆ lsm_fab_data

amrex::Vector<FabPtr> OceanSurf::lsm_fab_data
private

Referenced by Init(), and Lsm_Data_Ptr().

◆ LsmVarMap

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

Referenced by Init(), Lsm_Data_Ptr(), and Lsm_DataName().

◆ LsmVarName

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

Referenced by Init(), and Lsm_DataName().

◆ m_default_tsurf

amrex::Real OceanSurf::m_default_tsurf {amrex::Real(300.0)}
private

Referenced by Define(), and Init().

◆ m_dt

amrex::Real OceanSurf::m_dt {0.0}
private

Referenced by Advance(), and Init().

◆ m_lsm_geom

amrex::Geometry OceanSurf::m_lsm_geom
private

Referenced by Init(), and Lsm_Geom().

◆ m_lsm_size

int OceanSurf::m_lsm_size = 1
private

Referenced by Init(), and Lsm_Data_Size().

◆ m_ocean_state

const amrex::Vector<amrex::MultiFab*>* OceanSurf::m_ocean_state = nullptr
private

Referenced by OceanState(), and SetOceanState().


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