ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Utils.cpp File Reference
#include "ERF_Utils.H"
Include dependency graph for ERF_Utils.cpp:

Functions

void cons_to_prim (const MultiFab &cons_state, MultiFab &S_prim, int ng)
 
void make_qt (const MultiFab &cons_state, MultiFab &qt, int n_qstate_into_total)
 

Function Documentation

◆ cons_to_prim()

void cons_to_prim ( const MultiFab &  cons_state,
MultiFab &  S_prim,
int  ng 
)
7 {
8  BL_PROFILE("cons_to_prim()");
9 
10  int ncomp_prim = S_prim.nComp();
11 
12 #ifdef _OPENMP
13 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
14 #endif
15  for (MFIter mfi(cons_state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
16  {
17  const Box& gbx = mfi.growntilebox(ng);
18  const Array4<const Real>& cons_arr = cons_state.array(mfi);
19  const Array4< Real>& prim_arr = S_prim.array(mfi);
20 
21  //
22  // We may need > one ghost cells of prim in order to compute higher order advective terms
23  //
24  amrex::ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
25  {
26  Real rho = cons_arr(i,j,k,Rho_comp);
27  Real rho_theta = cons_arr(i,j,k,RhoTheta_comp);
28  prim_arr(i,j,k,PrimTheta_comp) = rho_theta / rho;
29  for (int n = 1; n < ncomp_prim; ++n) {
30  prim_arr(i,j,k,PrimTheta_comp + n) = cons_arr(i,j,k,RhoTheta_comp + n) / rho;
31  }
32  });
33  } // mfi
34 };
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:55
rho
Definition: ERF_InitCustomPert_Bubble.H:106
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ ng
Definition: ERF_Morrison.H:48

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

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

◆ make_qt()

void make_qt ( const MultiFab &  cons_state,
MultiFab &  qt,
int  n_qstate_into_total 
)
38 {
39  BL_PROFILE("make_qt()");
40 
41  // All moisture models are guaranteed to have RhoQ1_comp.
42  MultiFab::Copy(qt, cons_state, RhoQ1_comp, 0, 1, qt.nGrowVect());
43 
44  for (int n = 1; n < n_qstate_into_total; ++n) {
45  MultiFab::Add(qt, cons_state, RhoQ1_comp+n, 0, 1, qt.nGrowVect());
46  }
47 
48  MultiFab::Divide(qt, cons_state, Rho_comp, 0, 1, qt.nGrowVect());
49 }
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
@ qt
Definition: ERF_Kessler.H:28

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

Here is the caller graph for this function: