ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
derived Namespace Reference

Functions

void erf_derrhodivide (const Box &bx, FArrayBox &derfab, const FArrayBox &datfab, const int scalar_index)
 
void erf_dernull (const Box &, FArrayBox &, int, int, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dersoundspeed (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_dertemp (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_dertheta (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_derscalar (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_derKE (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_derQKE (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const Geometry &, Real, const int *, const int)
 
void erf_dervortx (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real, const int *, const int)
 
void erf_dervorty (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real, const int *, const int)
 
void erf_dervortz (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real, const int *, const int)
 
void erf_dermagvel (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &, amrex::Real, const int *, const int)
 
void erf_derrhodivide (const amrex::Box &bx, amrex::FArrayBox &derfab, const amrex::FArrayBox &datfab, const int scalar_index)
 
void erf_dernull (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dersoundspeed (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dertemp (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dertheta (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derscalar (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derKE (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derQKE (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 

Function Documentation

◆ erf_derKE() [1/2]

void derived::erf_derKE ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_derKE() [2/2]

void derived::erf_derKE ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define the kinetic energy KE by dividing (rho KE) by rho

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds KE @params[in] datfab array of data used to construct derived quantity

172 {
173  erf_derrhodivide(bx, derfab, datfab, RhoKE_comp);
174 }
#define RhoKE_comp
Definition: IndexDefines.H:14
void erf_derrhodivide(const Box &bx, FArrayBox &derfab, const FArrayBox &datfab, const int scalar_index)
Definition: Derive.cpp:18

Referenced by ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_dermagvel()

void derived::erf_dermagvel ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  ,
amrex::Real  ,
const int *  ,
const int   
)
292 {
293  AMREX_ALWAYS_ASSERT(dcomp == 0);
294  AMREX_ALWAYS_ASSERT(ncomp == 1);
295 
296  auto const dat = datfab.array(); // cell-centered velocity
297  auto tfab = derfab.array(); // cell-centered magvel
298 
299  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
300  {
301  Real u = dat(i,j,k,0);
302  Real v = dat(i,j,k,1);
303  Real w = dat(i,j,k,2);
304  tfab(i,j,k,dcomp) = std::sqrt(u*u + v*v + w*w);
305  });
306 }

Referenced by ERF::WritePlotFile().

Here is the caller graph for this function:

◆ erf_dernull() [1/2]

void derived::erf_dernull ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_dernull() [2/2]

void derived::erf_dernull ( const Box &  ,
FArrayBox &  ,
int  ,
int  ,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Placeholder function that does nothing

48 { }

Referenced by ERF::WritePlotFile().

Here is the caller graph for this function:

◆ erf_derQKE() [1/2]

void derived::erf_derQKE ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_derQKE() [2/2]

void derived::erf_derQKE ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define QKE by dividing (rho QKE) by rho

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds QKE @params[in] datfab array of data used to construct derived quantity

193 {
194  erf_derrhodivide(bx, derfab, datfab, RhoQKE_comp);
195 }
#define RhoQKE_comp
Definition: IndexDefines.H:15

Referenced by ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_derrhodivide() [1/2]

void derived::erf_derrhodivide ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
const amrex::FArrayBox &  datfab,
const int  scalar_index 
)

◆ erf_derrhodivide() [2/2]

void derived::erf_derrhodivide ( const Box &  bx,
FArrayBox &  derfab,
const FArrayBox &  datfab,
const int  scalar_index 
)

Function to define a derived quantity by dividing by density (analogous to cons_to_prim)

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity @params[in] datfab array of data used to construct derived quantity @params[in] scalar_index index of quantity to be divided by density

22 {
23  // This routine divides any cell-centered conserved quantity by density
24  auto const dat = datfab.array();
25  auto primitive = derfab.array();
26 
27  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
28  {
29  const Real rho = dat(i, j, k, Rho_comp);
30  const Real conserved = dat(i, j, k, scalar_index);
31  primitive(i,j,k) = conserved / rho;
32  });
33 }
#define Rho_comp
Definition: IndexDefines.H:12
@ rho
Definition: Kessler.H:30

Referenced by erf_derKE(), erf_derQKE(), erf_derscalar(), erf_dertheta(), and WriteBndryPlanes::write_planes().

Here is the caller graph for this function:

◆ erf_derscalar() [1/2]

void derived::erf_derscalar ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_derscalar() [2/2]

void derived::erf_derscalar ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define a scalar s by dividing (rho s) by rho

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds scalar s @params[in] datfab array of data used to construct derived quantity

151 {
152  erf_derrhodivide(bx, derfab, datfab, RhoScalar_comp);
153 }
#define RhoScalar_comp
Definition: IndexDefines.H:16

Referenced by ERF::ErrorEst(), and ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_dersoundspeed() [1/2]

void derived::erf_dersoundspeed ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_dersoundspeed() [2/2]

void derived::erf_dersoundspeed ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define the sound speed by calling an EOS routine

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds pressure @params[in] datfab array of data used to construct derived quantity

67 {
68  auto const dat = datfab.array();
69  auto cfab = derfab.array();
70 
71  // NOTE: we compute the soundspeed of dry air -- we do not account for any moisture effects here
72  Real qv = 0.;
73 
74  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
75  {
76  const Real rhotheta = dat(i, j, k, RhoTheta_comp);
77  const Real rho = dat(i, j, k, Rho_comp);
78  AMREX_ALWAYS_ASSERT(rhotheta > 0.);
79  cfab(i,j,k) = std::sqrt(Gamma * getPgivenRTh(rhotheta,qv) / rho);
80  });
81 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: EOS.H:62
constexpr amrex::Real Gamma
Definition: ERF_Constants.H:19
#define RhoTheta_comp
Definition: IndexDefines.H:13
@ qv
Definition: Kessler.H:36

Referenced by ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_dertemp() [1/2]

void derived::erf_dertemp ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_dertemp() [2/2]

void derived::erf_dertemp ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define the temperature by calling an EOS routine

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds pressure @params[in] datfab array of data used to construct derived quantity

100 {
101  auto const dat = datfab.array();
102  auto tfab = derfab.array();
103 
104  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
105  {
106  const Real rho = dat(i, j, k, Rho_comp);
107  const Real rhotheta = dat(i, j, k, RhoTheta_comp);
108  AMREX_ALWAYS_ASSERT(rhotheta > 0.);
109  tfab(i,j,k) = getTgivenRandRTh(rho,rhotheta);
110  });
111 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: EOS.H:29

Referenced by WriteBndryPlanes::write_planes(), and ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_dertheta() [1/2]

void derived::erf_dertheta ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  time,
const int *  bcrec,
const int  level 
)

◆ erf_dertheta() [2/2]

void derived::erf_dertheta ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Function to define the potential temperature by calling an EOS routine

@params[in] bx box on which to divide by density @params[out] derfab array of derived quantity – here it holds pressure @params[in] datfab array of data used to construct derived quantity

130 {
131  erf_derrhodivide(bx, derfab, datfab, RhoTheta_comp);
132 }

Referenced by ERF::ErrorEst(), and ERF::WritePlotFile().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ erf_dervortx()

void derived::erf_dervortx ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  ,
const int *  ,
const int   
)
208 {
209  AMREX_ALWAYS_ASSERT(dcomp == 0);
210  AMREX_ALWAYS_ASSERT(ncomp == 1);
211 
212  auto const dat = datfab.array(); // cell-centered velocity
213  auto tfab = derfab.array(); // cell-centered vorticity x-component
214 
215  const Real dy = geomdata.CellSize(1);
216  const Real dz = geomdata.CellSize(2);
217 
218  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
219  {
220  tfab(i,j,k,dcomp) = (dat(i,j+1,k,2) - dat(i,j-1,k,2)) / (2.0*dy) // dw/dy
221  - (dat(i,j,k+1,1) - dat(i,j,k-1,1)) / (2.0*dz); // dv/dz
222  });
223 }

Referenced by ERF::WritePlotFile().

Here is the caller graph for this function:

◆ erf_dervorty()

void derived::erf_dervorty ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  ,
const int *  ,
const int   
)
236 {
237  AMREX_ALWAYS_ASSERT(dcomp == 0);
238  AMREX_ALWAYS_ASSERT(ncomp == 1);
239 
240  auto const dat = datfab.array(); // cell-centered velocity
241  auto tfab = derfab.array(); // cell-centered vorticity y-component
242 
243  const Real dx = geomdata.CellSize(0);
244  const Real dz = geomdata.CellSize(2);
245 
246  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
247  {
248  tfab(i,j,k,dcomp) = (dat(i,j,k+1,0) - dat(i,j,k-1,0)) / (2.0*dz) // du/dz
249  - (dat(i+1,j,k,2) - dat(i-1,j,k,2)) / (2.0*dx); // dw/dx
250  });
251 }

Referenced by ERF::WritePlotFile().

Here is the caller graph for this function:

◆ erf_dervortz()

void derived::erf_dervortz ( const amrex::Box &  bx,
amrex::FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const amrex::FArrayBox &  datfab,
const amrex::Geometry &  geomdata,
amrex::Real  ,
const int *  ,
const int   
)
264 {
265  AMREX_ALWAYS_ASSERT(dcomp == 0);
266  AMREX_ALWAYS_ASSERT(ncomp == 1);
267 
268  auto const dat = datfab.array(); // cell-centered velocity
269  auto tfab = derfab.array(); // cell-centered vorticity z-component
270 
271  const Real dx = geomdata.CellSize(0);
272  const Real dy = geomdata.CellSize(2);
273 
274  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
275  {
276  tfab(i,j,k,dcomp) = (dat(i+1,j,k,1) - dat(i-1,j,k,1)) / (2.0*dx) // dv/dx
277  - (dat(i,j+1,k,0) - dat(i,j-1,k,0)) / (2.0*dy); // du/dy
278  });
279 }

Referenced by ERF::WritePlotFile().

Here is the caller graph for this function: