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 FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dersoundspeed (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dertemp (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dermoisttemp (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dertheta (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_derscalar (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_derKE (const Box &bx, FArrayBox &derfab, int, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dervortx (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &geomdata, Real, const int *, const int)
 
void erf_dervorty (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &geomdata, Real, const int *, const int)
 
void erf_dervortz (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &, const Geometry &geomdata, Real, const int *, const int)
 
void erf_derenstrophysq (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &geomdata, Real, const int *, const int)
 
void erf_dermagvel (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dermagvelsq (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_derreflectivity (const Box &bx, FArrayBox &derfab, int dcomp, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_dermaxreflectivity (const Box &bx, FArrayBox &derfab, int dcomp, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &, Real, const int *, const int)
 
void erf_derlocalhelicity (const Box &bx, FArrayBox &derfab, int dcomp, int, const FArrayBox &datfab, const FArrayBox &, const Geometry &geomdata, Real, const int *, const int)
 
void erf_derhelicity (const Box &bx, FArrayBox &derfab, int dcomp, int, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &geomdata, Real, const int *, const int)
 
void erf_derprecipitable (const Box &bx, FArrayBox &derfab, int dcomp, int, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &, Real, const int *, const int)
 
void erf_dermucape (const Box &bx, FArrayBox &derfab, int dcomp, int ncomp, const FArrayBox &datfab, const FArrayBox &zcc_fab, const Geometry &, 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::FArrayBox &zfab, 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::FArrayBox &zfab, 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::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dermoisttemp (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, 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::FArrayBox &zfab, 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::FArrayBox &zfab, 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::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dervortx (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dervorty (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dervortz (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derenstrophysq (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dermagvel (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dermagvelsq (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derreflectivity (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dermaxreflectivity (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derhelicity (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derlocalhelicity (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_derprecipitable (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 
void erf_dermucape (const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::FArrayBox &zfab, const amrex::Geometry &geomdata, amrex::Real time, const int *bcrec, const int level)
 

Function Documentation

◆ erf_derenstrophysq() [1/2]

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

◆ erf_derenstrophysq() [2/2]

void derived::erf_derenstrophysq ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
453 {
454  AMREX_ALWAYS_ASSERT(dcomp == 0);
455  AMREX_ALWAYS_ASSERT(ncomp == 1);
456 
457  auto const dat = datfab.array(); // cell-centered velocity
458  auto tfab = derfab.array(); // cell-centered vorticity x-component
459  auto z_arr = zcc_fab.array(); // cell-centered height z
460 
461  const Real two_dx = two * geomdata.CellSize(0);
462  const Real two_dy = two * geomdata.CellSize(1);
463 
464  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
465  {
466  Real two_dz = z_arr(i,j,k+1) - z_arr(i,j,k-1);
467 
468  Real vortx = (dat(i,j+1,k,2) - dat(i,j-1,k,2)) / two_dy // dw/dy
469  -(dat(i,j,k+1,1) - dat(i,j,k-1,1)) / two_dz; // dv/dz
470  Real vorty = (dat(i,j,k+1,0) - dat(i,j,k-1,0)) / two_dz // du/dz
471  -(dat(i+1,j,k,2) - dat(i-1,j,k,2)) / two_dx; // dw/dx
472  Real vortz = (dat(i+1,j,k,1) - dat(i-1,j,k,1)) / two_dx // dv/dx
473  -(dat(i,j+1,k,0) - dat(i,j-1,k,0)) / two_dy; // du/dy
474 
475  tfab(i,j,k,dcomp) = vortx*vortx + vorty*vorty + vortz*vortz;
476  });
477 }
constexpr amrex::Real two
Definition: ERF_Constants.H:8
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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);} } })
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by ERF::sum_derived_quantities().

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

◆ erf_derhelicity() [1/2]

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

◆ erf_derhelicity() [2/2]

void derived::erf_derhelicity ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
656 {
657  AMREX_ALWAYS_ASSERT(dcomp == 0);
658 
659  auto const dat = datfab.array(); // cell-centered velocity
660  auto dfab = derfab.array(); // integral of local helicity
661  auto z_arr = zcc_fab.array(); // cell-centered height z
662 
663  const Real dx = geomdata.CellSize(0);
664  const Real dy = geomdata.CellSize(1);
665 
666  // Collapse to i,j box (ignore vertical for now)
667  Box b2d = bx;
668  b2d.setSmall(2,0);
669  b2d.setBig(2,0);
670 
671  ParallelFor(b2d, [=] AMREX_GPU_DEVICE(int i, int j, int ) noexcept
672  {
673  Real int_hel = Real(0.0);
674  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k)
675  {
676  Real z = z_arr(i,j,k);
677 
678  // Helicity is defined as integral from 2km to 5km in vertical
679  if (z > Real(2000.0) && z < Real(5000.0)) {
680 
681  Real z_hi = Real(0.5) * (z_arr(i,j,k) + z_arr(i,j,k+1));
682  Real z_lo = Real(0.5) * (z_arr(i,j,k) + z_arr(i,j,k-1));
683  Real dz = z_hi - z_lo;
684 
685  Real vortz = (dat(i+1,j,k,1) - dat(i-1,j,k,1)) / (2.0*dx) // dv/dx
686  - (dat(i,j+1,k,0) - dat(i,j-1,k,0)) / (2.0*dy); // du/dy
687  Real w = dat(i,j,k,2); // vertical velocity
688 
689  int_hel += vortz * w * dz;
690  }
691  }
692 
693  // Store vertical integral into *all* levels for this (i,j)
694  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k) {
695  dfab(i, j, k, dcomp) = int_hel;
696  }
697  });
698 }
const Real dy
Definition: ERF_InitCustomPert_ABL.H:24
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25

Referenced by ERF::Write3DPlotFile().

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

◆ 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::FArrayBox &  zfab,
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 FArrayBox &  ,
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

352 {
353  erf_derrhodivide(bx, derfab, datfab, RhoKE_comp);
354 }
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
void erf_derrhodivide(const Box &bx, FArrayBox &derfab, const FArrayBox &datfab, const int scalar_index)
Definition: ERF_Derive.cpp:168

Referenced by ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

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

◆ erf_derlocalhelicity() [1/2]

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

◆ erf_derlocalhelicity() [2/2]

void derived::erf_derlocalhelicity ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
625 {
626  AMREX_ALWAYS_ASSERT(dcomp == 0);
627 
628  auto const dat = datfab.array(); // cell-centered velocity
629  auto dfab = derfab.array(); // cell-centered local helicity
630 
631  const Real two_dx = Real(2.0)*geomdata.CellSize(0);
632  const Real two_dy = Real(2.0)*geomdata.CellSize(1);
633 
634  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
635  {
636  Real vortz = (dat(i+1,j,k,1) - dat(i-1,j,k,1)) / two_dx // dv/dx
637  - (dat(i,j+1,k,0) - dat(i,j-1,k,0)) / two_dy; // du/dy
638  Real w = dat(i,j,k,2);
639 
640  // Helicity
641  dfab(i,j,k,dcomp) = vortz * w;
642  });
643 }

Referenced by ERF::Write3DPlotFile().

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

◆ erf_dermagvel() [1/2]

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

◆ erf_dermagvel() [2/2]

void derived::erf_dermagvel ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
490 {
491  AMREX_ALWAYS_ASSERT(dcomp == 0);
492  AMREX_ALWAYS_ASSERT(ncomp == 1);
493 
494  auto const dat = datfab.array(); // cell-centered velocity
495  auto tfab = derfab.array(); // cell-centered magvel
496 
497  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
498  {
499  Real u = dat(i,j,k,0);
500  Real v = dat(i,j,k,1);
501  Real w = dat(i,j,k,2);
502  tfab(i,j,k,dcomp) = std::sqrt(u*u + v*v + w*w);
503  });
504 }

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

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

◆ erf_dermagvelsq() [1/2]

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

◆ erf_dermagvelsq() [2/2]

void derived::erf_dermagvelsq ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
517 {
518  AMREX_ALWAYS_ASSERT(dcomp == 0);
519  AMREX_ALWAYS_ASSERT(ncomp == 1);
520 
521  auto const dat = datfab.array(); // cell-centered velocity
522  auto tfab = derfab.array(); // cell-centered magvel
523 
524  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
525  {
526  Real u = dat(i,j,k,0);
527  Real v = dat(i,j,k,1);
528  Real w = dat(i,j,k,2);
529  tfab(i,j,k,dcomp) = u*u + v*v + w*w;
530  });
531 }

Referenced by ERF::sum_derived_quantities().

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

◆ erf_dermaxreflectivity() [1/2]

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

◆ erf_dermaxreflectivity() [2/2]

void derived::erf_dermaxreflectivity ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
576 {
577  AMREX_ALWAYS_ASSERT(dcomp == 0);
578 
579  auto const dat = datfab.array(); // cell-centered state vector
580  auto rfab = derfab.array(); // cell-centered max reflectivity
581 
582  // Collapse to i,j box (ignore vertical for now)
583  Box b2d = bx;
584  b2d.setSmall(2,0);
585  b2d.setBig(2,0);
586 
587  ParallelFor(b2d, [=] AMREX_GPU_DEVICE(int i, int j, int) noexcept {
588 
589  Real max_dbz = -1.0e30;
590 
591  // find max reflectivity over k
592  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k) {
593 
594  Real rho = dat(i,j,k,Rho_comp);
595  Real qv = std::max(Real(0.0),dat(i,j,k,RhoQ1_comp)/rho);
596  Real qpr = std::max(Real(0.0),dat(i,j,k,RhoQ4_comp)/rho);
597  Real qps = std::max(Real(0.0),dat(i,j,k,RhoQ5_comp)/rho);
598  Real qpg = std::max(Real(0.0),dat(i,j,k,RhoQ6_comp)/rho);
599 
600  Real temp = getTgivenRandRTh(rho, dat(i,j,k,RhoTheta_comp), qv);
601 
602  Real dbz = compute_max_reflectivity_dbz(rho, temp, qpr, qps, qpg,
603  1, 1, 1, 1);
604  max_dbz = amrex::max(max_dbz, dbz);
605  }
606 
607  // Store max_dbz into *all* levels for this (i,j)
608  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k) {
609  rfab(i, j, k, dcomp) = max_dbz;
610  }
611  });
612 }
for(int i=0;i< m_num_species;i++)
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:48

Referenced by ERF::Write3DPlotFile().

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

◆ erf_dermoisttemp() [1/2]

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

◆ erf_dermoisttemp() [2/2]

void derived::erf_dermoisttemp ( const Box &  bx,
FArrayBox &  derfab,
int  ,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
276 {
277  auto const dat = datfab.array();
278  auto tfab = derfab.array();
279 
280  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
281  {
282  const Real rho = dat(i, j, k, Rho_comp);
283  const Real rhotheta = dat(i, j, k, RhoTheta_comp);
284  AMREX_ALWAYS_ASSERT(rhotheta > Real(0.0));
285  const Real qv = dat(i, j, k, RhoQ1_comp) / rho;
286  tfab(i,j,k) = getTgivenRandRTh(rho,rhotheta,qv);
287  });
288 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
rho
Definition: ERF_InitCustomPert_Bubble.H:106
@ qv
Definition: ERF_Kessler.H:28

Referenced by ERF::Write3DPlotFile(), WriteBndryPlanes::write_planes(), and ERF::WriteSubvolume().

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

◆ erf_dermucape() [1/2]

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

◆ erf_dermucape() [2/2]

void derived::erf_dermucape ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
758 {
759  AMREX_ALWAYS_ASSERT(dcomp == 0);
760  AMREX_ALWAYS_ASSERT(ncomp == 1);
761 
762  auto const dat = datfab.array();
763  auto dfab = derfab.array();
764  auto z_arr = zcc_fab.array();
765  const int ncons = datfab.nComp();
766 
767  Box b2d = bx;
768  b2d.setSmall(2,0);
769  b2d.setBig(2,0);
770 
771  ParallelFor(b2d, [=] AMREX_GPU_DEVICE(int i, int j, int) noexcept
772  {
773  Real mucape = Real(0);
774  int klo = bx.smallEnd(2);
775  int khi = bx.bigEnd(2);
776 
777  if (ncons > RhoQ1_comp) {
778  Real p_sfc = Real(0);
779  for (int k = klo; k <= khi; ++k) {
780  Real rho = dat(i,j,k,Rho_comp);
781  if (rho <= Real(0)) { continue; }
782  Real qv = amrex::max(Real(0), dat(i,j,k,RhoQ1_comp) / rho);
783  p_sfc = amrex::max(p_sfc, getPgivenRTh(dat(i,j,k,RhoTheta_comp), qv));
784  }
785 
786  if (p_sfc > mucape_search_depth_pa()) {
787  Real p_search_min = p_sfc - mucape_search_depth_pa();
788 
789  for (int ks = klo; ks < khi; ++ks) {
790  Real rho_src = dat(i,j,ks,Rho_comp);
791  if (rho_src <= Real(0)) { continue; }
792 
793  Real qv_src = amrex::max(Real(0), dat(i,j,ks,RhoQ1_comp) / rho_src);
794  Real rt_src = dat(i,j,ks,RhoTheta_comp);
795  Real p_src = getPgivenRTh(rt_src, qv_src);
796 
797  if (p_src < p_search_min) { continue; }
798 
799  Real T_src = getTgivenRandRTh(rho_src, rt_src, qv_src);
800  if (T_src <= Real(0)) { continue; }
801 
802  Real Td_src = mucape_dewpoint_temperature(p_src, qv_src, T_src);
803  Real Tlcl = mucape_lcl_temperature(T_src, Td_src);
804  Real plcl = p_src * std::pow(Tlcl / T_src, Cp_d / R_d);
805  plcl = amrex::min(p_src, amrex::max(plcl, mucape_min_pressure_pa()));
806 
807  Real theta_src = getThgivenTandP(T_src, p_src, R_d / Cp_d);
808 
809  Real candidate_cape = Real(0);
810  Real z_prev = z_arr(i,j,ks);
811  Real b_prev = Real(0);
812 
813  bool saturated = false;
814  Real T_sat_prev = Tlcl;
815  Real p_sat_prev = plcl;
816 
817  for (int k = ks + 1; k <= khi; ++k) {
818  Real rho_env = dat(i,j,k,Rho_comp);
819  if (rho_env <= Real(0)) { continue; }
820 
821  Real qv_env = amrex::max(Real(0), dat(i,j,k,RhoQ1_comp) / rho_env);
822  Real rt_env = dat(i,j,k,RhoTheta_comp);
823  Real p_env = getPgivenRTh(rt_env, qv_env);
824  Real T_env = getTgivenRandRTh(rho_env, rt_env, qv_env);
825  Real Tv_env = mucape_virtual_temperature(T_env, qv_env);
826  Real z_env = z_arr(i,j,k);
827 
828  Real T_parcel;
829  Real qv_parcel;
830 
831  if (p_env >= plcl) {
832  T_parcel = getTgivenPandTh(p_env, theta_src, R_d / Cp_d);
833  qv_parcel = qv_src;
834  } else {
835  if (!saturated) {
836  saturated = true;
837  T_sat_prev = Tlcl;
838  p_sat_prev = plcl;
839  }
840 
841  if (p_env < p_sat_prev) {
842  T_sat_prev = mucape_integrate_saturated_temperature(T_sat_prev, p_sat_prev, p_env);
843  p_sat_prev = p_env;
844  }
845 
846  T_parcel = T_sat_prev;
847  qv_parcel = mucape_qsat(T_parcel, p_env);
848  }
849 
850  Real Tv_parcel = mucape_virtual_temperature(T_parcel, qv_parcel);
851  Real buoyancy = CONST_GRAV * (Tv_parcel - Tv_env) /
852  amrex::max(Tv_env, mucape_min_temperature());
853 
854  candidate_cape += mucape_positive_area(b_prev, buoyancy, z_env - z_prev);
855  b_prev = buoyancy;
856  z_prev = z_env;
857  }
858 
859  mucape = amrex::max(mucape, candidate_cape);
860  }
861  }
862  }
863 
864  for (int k = klo; k <= khi; ++k) {
865  dfab(i,j,k,dcomp) = mucape;
866  }
867  });
868 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:22
constexpr amrex::Real CONST_GRAV
Definition: ERF_Constants.H:31
constexpr amrex::Real R_d
Definition: ERF_Constants.H:20
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenPandTh(const amrex::Real P, const amrex::Real th, const amrex::Real rdOcp)
Definition: ERF_EOS.H:32
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21

Referenced by ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

Here is the call graph for this function:
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::FArrayBox &  zfab,
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 FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)

Placeholder function that does nothing

199 { }

Referenced by ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

Here is the caller graph for this function:

◆ erf_derprecipitable() [1/2]

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

◆ erf_derprecipitable() [2/2]

void derived::erf_derprecipitable ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
711 {
712  AMREX_ALWAYS_ASSERT(dcomp == 0);
713 
714  auto const dat = datfab.array(); // cell-centered state vector
715  auto dfab = derfab.array(); // integral of qv to define precipitable water
716 
717  // Collapse to i,j box (ignore vertical for now)
718  Box b2d = bx;
719  b2d.setSmall(2,0);
720  b2d.setBig(2,0);
721 
722  auto z_arr = zcc_fab.array(); // cell-centered height z
723 
724  ParallelFor(b2d, [=] AMREX_GPU_DEVICE(int i, int j, int) noexcept
725  {
726 
727  Real int_qv = Real(0.0);
728 
729  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k)
730  {
731  Real z_hi = Real(0.5) * (z_arr(i,j,k) + z_arr(i,j,k+1));
732  Real z_lo = Real(0.5) * (z_arr(i,j,k) + z_arr(i,j,k-1));
733  Real dz = z_hi - z_lo;
734 
735  Real rhoQ1 = dat(i, j, k, RhoQ1_comp);
736 
737  int_qv += rhoQ1 * dz;
738  }
739 
740  // Store vertical integral into *all* levels for this (i,j)
741  for (int k = bx.smallEnd(2); k <= bx.bigEnd(2); ++k) {
742  dfab(i, j, k, dcomp) = int_qv;
743  }
744  });
745 }

Referenced by ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

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

◆ erf_derreflectivity() [1/2]

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

◆ erf_derreflectivity() [2/2]

void derived::erf_derreflectivity ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  ,
Real  ,
const int *  ,
const int   
)
544 {
545  AMREX_ALWAYS_ASSERT(dcomp == 0);
546 
547  auto const dat = datfab.array(); // cell-centered state vector
548  auto rfab = derfab.array(); // cell-centered reflectivity
549 
550  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
551  {
552  Real rho = dat(i,j,k,Rho_comp);
553  Real qv = std::max(Real(0.0),dat(i,j,k,RhoQ1_comp)/rho);
554  Real qpr = std::max(Real(0.0),dat(i,j,k,RhoQ4_comp)/rho);
555  Real qps = std::max(Real(0.0),dat(i,j,k,RhoQ5_comp)/rho);
556  Real qpg = std::max(Real(0.0),dat(i,j,k,RhoQ6_comp)/rho);
557 
558  Real temp = getTgivenRandRTh(rho, dat(i,j,k,RhoTheta_comp), qv);
559 
560  rfab(i, j, k, dcomp) = compute_max_reflectivity_dbz(rho, temp, qpr, qps, qpg,
561  1, 1, 1, 1);
562  });
563 }
#define RhoQ4_comp
Definition: ERF_IndexDefines.H:45
#define RhoQ6_comp
Definition: ERF_IndexDefines.H:47
#define RhoQ5_comp
Definition: ERF_IndexDefines.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real compute_max_reflectivity_dbz(amrex::Real rho_air, amrex::Real tmk, amrex::Real qra, amrex::Real qsn, amrex::Real qgr, int in0r, int in0s, int in0g, int iliqskin)
Definition: ERF_StormDiagnostics.H:13
@ qpg
Definition: ERF_Morrison.H:41
@ qps
Definition: ERF_Morrison.H:40
@ qpr
Definition: ERF_Morrison.H:39

Referenced by ERF::Write3DPlotFile().

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

172 {
173  // This routine divides any cell-centered conserved quantity by density
174  auto const dat = datfab.array();
175  auto primitive = derfab.array();
176 
177  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
178  {
179  const Real rho = dat(i, j, k, Rho_comp);
180  const Real conserved = dat(i, j, k, scalar_index);
181  primitive(i,j,k) = conserved / rho;
182  });
183 }

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

Here is the call graph for this function:
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::FArrayBox &  zfab,
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 FArrayBox &  ,
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

330 {
331  erf_derrhodivide(bx, derfab, datfab, RhoScalar_comp);
332 }
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40

Referenced by ERF::ErrorEst(), ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

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::FArrayBox &  zfab,
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 FArrayBox &  ,
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

219 {
220  auto const dat = datfab.array();
221  auto cfab = derfab.array();
222 
223  // NOTE: we compute the soundspeed of dry air -- we do not account for any moisture effects here
224  Real qv = Real(0.0);
225 
226  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
227  {
228  const Real rhotheta = dat(i, j, k, RhoTheta_comp);
229  const Real rho = dat(i, j, k, Rho_comp);
230  AMREX_ALWAYS_ASSERT(rhotheta > 0);
231  cfab(i,j,k) = std::sqrt(Gamma * getPgivenRTh(rhotheta,qv) / rho);
232  });
233 }
constexpr amrex::Real Gamma
Definition: ERF_Constants.H:29

Referenced by ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

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::FArrayBox &  zfab,
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 FArrayBox &  ,
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

253 {
254  auto const dat = datfab.array();
255  auto tfab = derfab.array();
256 
257  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
258  {
259  const Real rho = dat(i, j, k, Rho_comp);
260  const Real rhotheta = dat(i, j, k, RhoTheta_comp);
261  AMREX_ALWAYS_ASSERT(rhotheta > Real(0.0));
262  tfab(i,j,k) = getTgivenRandRTh(rho,rhotheta);
263  });
264 }

Referenced by ERF::Write3DPlotFile(), WriteBndryPlanes::write_planes(), and ERF::WriteSubvolume().

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::FArrayBox &  zfab,
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 FArrayBox &  ,
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

308 {
309  erf_derrhodivide(bx, derfab, datfab, RhoTheta_comp);
310 }

Referenced by ERF::ErrorEst(), ERF::Write3DPlotFile(), and ERF::WriteSubvolume().

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

◆ erf_dervortx() [1/2]

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

◆ erf_dervortx() [2/2]

void derived::erf_dervortx ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
367 {
368  AMREX_ALWAYS_ASSERT(dcomp == 0);
369  AMREX_ALWAYS_ASSERT(ncomp == 1);
370 
371  auto const dat = datfab.array(); // cell-centered velocity
372  auto tfab = derfab.array(); // cell-centered vorticity x-component
373  auto z_arr = zcc_fab.array(); // cell-centered height z
374 
375  const Real two_dy = two * geomdata.CellSize(1);
376 
377  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
378  {
379  Real two_dz = z_arr(i,j,k+1) - z_arr(i,j,k-1);
380  tfab(i,j,k,dcomp) = (dat(i,j+1,k,2) - dat(i,j-1,k,2)) / two_dy // dw/dy
381  - (dat(i,j,k+1,1) - dat(i,j,k-1,1)) / two_dz; // dv/dz
382  });
383 }

Referenced by ERF::Write3DPlotFile().

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

◆ erf_dervorty() [1/2]

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

◆ erf_dervorty() [2/2]

void derived::erf_dervorty ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  zcc_fab,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
396 {
397  AMREX_ALWAYS_ASSERT(dcomp == 0);
398  AMREX_ALWAYS_ASSERT(ncomp == 1);
399 
400  auto const dat = datfab.array(); // cell-centered velocity
401  auto tfab = derfab.array(); // cell-centered vorticity y-component
402  auto z_arr = zcc_fab.array(); // cell-centered height z
403 
404  const Real two_dx = two * geomdata.CellSize(0);
405 
406  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
407  {
408  Real two_dz = z_arr(i,j,k+1) - z_arr(i,j,k-1);
409  tfab(i,j,k,dcomp) = (dat(i,j,k+1,0) - dat(i,j,k-1,0)) / two_dz // du/dz
410  - (dat(i+1,j,k,2) - dat(i-1,j,k,2)) / two_dx; // dw/dx
411  });
412 }

Referenced by ERF::Write3DPlotFile().

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

◆ erf_dervortz() [1/2]

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

◆ erf_dervortz() [2/2]

void derived::erf_dervortz ( const Box &  bx,
FArrayBox &  derfab,
int  dcomp,
int  ncomp,
const FArrayBox &  datfab,
const FArrayBox &  ,
const Geometry &  geomdata,
Real  ,
const int *  ,
const int   
)
425 {
426  AMREX_ALWAYS_ASSERT(dcomp == 0);
427  AMREX_ALWAYS_ASSERT(ncomp == 1);
428 
429  auto const dat = datfab.array(); // cell-centered velocity
430  auto tfab = derfab.array(); // cell-centered vorticity z-component
431 
432  const Real dx = geomdata.CellSize(0);
433  const Real dy = geomdata.CellSize(2);
434 
435  ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
436  {
437  tfab(i,j,k,dcomp) = (dat(i+1,j,k,1) - dat(i-1,j,k,1)) / (two*dx) // dv/dx
438  - (dat(i,j+1,k,0) - dat(i,j-1,k,0)) / (two*dy); // du/dy
439  });
440 }

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

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