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

#include <ERF_ForestDrag.H>

Collaboration diagram for ForestDrag:

Public Member Functions

 ForestDrag (std::string forestfile)
 
 ~ForestDrag ()=default
 
void define_drag_field (const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, amrex::Geometry &geom, amrex::MultiFab *z_phys_cc, amrex::MultiFab *z_phys_nd)
 
amrex::MultiFab * get_drag_field ()
 

Private Attributes

amrex::Vector< amrex::Realm_type_forest
 
amrex::Vector< amrex::Realm_x_forest
 
amrex::Vector< amrex::Realm_y_forest
 
amrex::Vector< amrex::Realm_height_forest
 
amrex::Vector< amrex::Realm_diameter_forest
 
amrex::Vector< amrex::Realm_cd_forest
 
amrex::Vector< amrex::Realm_lai_forest
 
amrex::Vector< amrex::Realm_laimax_forest
 
std::unique_ptr< amrex::MultiFab > m_forest_drag
 

Constructor & Destructor Documentation

◆ ForestDrag()

ForestDrag::ForestDrag ( std::string  forestfile)
explicit
11 {
12  std::ifstream file(forestfile, std::ios::in);
13  if (!file.good()) {
14  Abort("Cannot find forest file: " + forestfile);
15  }
16  // TreeType xc yc height diameter cd lai laimax
17  Real value1, value2, value3, value4, value5, value6, value7, value8;
18  while (file >> value1 >> value2 >> value3 >> value4 >> value5 >> value6 >>
19  value7 >> value8) {
20  m_type_forest.push_back(value1);
21  m_x_forest.push_back(value2);
22  m_y_forest.push_back(value3);
23  m_height_forest.push_back(value4);
24  m_diameter_forest.push_back(value5);
25  m_cd_forest.push_back(value6);
26  m_lai_forest.push_back(value7);
27  m_laimax_forest.push_back(value8);
28  }
29  file.close();
30 }
struct @22 in
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::Vector< amrex::Real > m_diameter_forest
Definition: ERF_ForestDrag.H:35
amrex::Vector< amrex::Real > m_laimax_forest
Definition: ERF_ForestDrag.H:38
amrex::Vector< amrex::Real > m_lai_forest
Definition: ERF_ForestDrag.H:37
amrex::Vector< amrex::Real > m_x_forest
Definition: ERF_ForestDrag.H:32
amrex::Vector< amrex::Real > m_height_forest
Definition: ERF_ForestDrag.H:34
amrex::Vector< amrex::Real > m_cd_forest
Definition: ERF_ForestDrag.H:36
amrex::Vector< amrex::Real > m_type_forest
Definition: ERF_ForestDrag.H:31
amrex::Vector< amrex::Real > m_y_forest
Definition: ERF_ForestDrag.H:33

◆ ~ForestDrag()

ForestDrag::~ForestDrag ( )
default

Member Function Documentation

◆ define_drag_field()

void ForestDrag::define_drag_field ( const amrex::BoxArray &  ba,
const amrex::DistributionMapping &  dm,
amrex::Geometry &  geom,
amrex::MultiFab *  z_phys_cc,
amrex::MultiFab *  z_phys_nd 
)
38 {
39  // Geometry params
40  const auto& dx = geom.CellSizeArray();
41  const auto& prob_lo = geom.ProbLoArray();
42 
43  bool all_boxes_touch_bottom = true;
44  for (int i = 0; i < ba.size(); i++) {
45  if (ba[i].smallEnd(2) != geom.ProbLo(2)) {
46  all_boxes_touch_bottom = false;
47  }
48  }
49  AMREX_ALWAYS_ASSERT(all_boxes_touch_bottom);
50 
51  // Allocate the forest drag MF
52  // NOTE: 1 ghost cell for averaging to faces
53  m_forest_drag.reset();
54  m_forest_drag = std::make_unique<MultiFab>(ba,dm,1,1);
55  m_forest_drag->setVal(0.);
56 
57  // Loop over forest types and pre-compute factors
58  for (unsigned ii = 0; ii < m_x_forest.size(); ++ii)
59  {
60  // Expose CPU data for GPU capture
61  Real af; // Depends upon the type of forest (tf)
62  Real treeZm = zero; // Only for forest type 2
63  int tf = int(m_type_forest[ii]);
64  Real hf = m_height_forest[ii];
65  Real xf = m_x_forest[ii];
66  Real yf = m_y_forest[ii];
67  Real df = m_diameter_forest[ii];
68  Real cdf = m_cd_forest[ii];
69  Real laif = m_lai_forest[ii];
70  Real laimaxf = m_laimax_forest[ii];
71  if (tf == 1) {
72  // Constant factor
73  af = laif / hf;
74  } else {
75  // Discretize integral with 100 points and pre-compute
76  int nk = 100;
77  Real ztree = 0;
78  Real expFun = 0;
79  Real ratio = 0;
80  const Real dz = hf / Real(nk);
81  treeZm = laimaxf * hf;
82  for (int k(0); k<nk; ++k) {
83  ratio = (hf - treeZm) / (hf - ztree);
84  if (ztree < treeZm) {
85  expFun += std::pow(ratio, Real(6.0)) *
86  std::exp(6 * (1 - ratio));
87  } else {
88  expFun += std::pow(ratio, myhalf) *
89  std::exp(myhalf * (1 - ratio));
90  }
91  ztree += dz;
92  }
93  af = laif / (expFun * dz);
94  }
95 
96  // Set the forest drag data
97  for (MFIter mfi(*m_forest_drag); mfi.isValid(); ++mfi) {
98  Box gtbx = mfi.growntilebox();
99  const Array4<Real>& levelDrag = m_forest_drag->array(mfi);
100  const Array4<const Real>& z_cc = z_phys_cc->const_array(mfi);
101  const Array4<const Real>& z_nd = z_phys_nd->const_array(mfi);
102 
103  ParallelFor(gtbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
104  {
105  // Physical positions of cell-centers
106  const Real x = prob_lo[0] + (i + myhalf) * dx[0];
107  const Real y = prob_lo[1] + (j + myhalf) * dx[1];
108 
109  // "z" is measured as distance from cell center to ground
110  const Real z_sfc = fourth * ( z_nd(i,j ,0) + z_nd(i+1,j ,0)
111  +z_nd(i,j+1,0) + z_nd(i+1,j+1,0));
112  const Real z = std::max((z_cc(i,j,k)-z_sfc),amrex::Real(0));
113 
114  // Proximity to the forest
115  const Real radius = std::sqrt((x - xf) * (x - xf) +
116  (y - yf) * (y - yf));
117 
118  // Hit for canopy region
119  Real factor = 1;
120  if ((z <= hf) && (radius <= (myhalf * df))) {
121  if (tf == 2) {
122  Real ratio = (hf - treeZm) / (hf - z);
123  if (z < treeZm) {
124  factor = std::pow(ratio, Real(6.0)) *
125  std::exp(Real(6.0) * (one - ratio));
126  } else if (z <= hf) {
127  factor = std::pow(ratio, myhalf) *
128  std::exp(myhalf * (one - ratio));
129  }
130  }
131  levelDrag(i, j, k) = cdf * af * factor;
132  }
133  });
134  } // mfi
135  } // ii (forest type)
136 
137  // Fillboundary for periodic ghost cell copy
138  m_forest_drag->FillBoundary(geom.periodicity());
139 
140 } // init_drag_field
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
const amrex::Real * prob_lo
Definition: ERF_InitCustomPert_IsentropicVortex.H:16
std::unique_ptr< amrex::MultiFab > m_forest_drag
Definition: ERF_ForestDrag.H:39
Here is the call graph for this function:

◆ get_drag_field()

amrex::MultiFab* ForestDrag::get_drag_field ( )
inline
28 { return m_forest_drag.get(); }

Member Data Documentation

◆ m_cd_forest

amrex::Vector<amrex::Real> ForestDrag::m_cd_forest
private

◆ m_diameter_forest

amrex::Vector<amrex::Real> ForestDrag::m_diameter_forest
private

◆ m_forest_drag

std::unique_ptr<amrex::MultiFab> ForestDrag::m_forest_drag
private

Referenced by get_drag_field().

◆ m_height_forest

amrex::Vector<amrex::Real> ForestDrag::m_height_forest
private

◆ m_lai_forest

amrex::Vector<amrex::Real> ForestDrag::m_lai_forest
private

◆ m_laimax_forest

amrex::Vector<amrex::Real> ForestDrag::m_laimax_forest
private

◆ m_type_forest

amrex::Vector<amrex::Real> ForestDrag::m_type_forest
private

◆ m_x_forest

amrex::Vector<amrex::Real> ForestDrag::m_x_forest
private

◆ m_y_forest

amrex::Vector<amrex::Real> ForestDrag::m_y_forest
private

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