ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TI_slow_rhs_post.H
Go to the documentation of this file.
1 #include "ERF_SrcHeaders.H"
2 
3  auto slow_rhs_fun_post = [&](Vector<MultiFab>& S_rhs,
4  Vector<MultiFab>& S_old,
5  Vector<MultiFab>& S_new,
6  Vector<MultiFab>& S_data,
7  const Real old_step_time,
8  const Real old_stage_time,
9  const Real new_stage_time,
10  const int nrk)
11  {
12  // Note that the "old" and "new" metric terms correspond to
13  // t^n and the RK stage (either t^*, t^** or t^{n+1} that this source
14  // will be used to advance to
15  Real slow_dt = new_stage_time - old_step_time;
16 
17  if (verbose) amrex::Print() << "Time integration of scalars at level " << level
18  << std::setprecision(timeprecision)
19  << " from " << old_step_time << " to " << new_stage_time
20  << " with dt = " << slow_dt
21  << " using RHS created at " << old_stage_time << std::endl;
22 
23 #if defined(ERF_USE_NETCDF)
24  bool moist_set_rhs = false;
25  if ( solverChoice.use_real_bcs && (level==0) &&
26  (solverChoice.moisture_type != MoistureType::None) )
27  {
28  moist_set_rhs = true;
29  }
30 #endif
31 
32  const GpuArray<Real, AMREX_SPACEDIM> dxInv = fine_geom.InvCellSizeArray();
33 
34  // *************************************************************************
35  // Set up flux registers if using two_way coupling
36  // *************************************************************************
37  YAFluxRegister* fr_as_crse = nullptr;
38  YAFluxRegister* fr_as_fine = nullptr;
39  if (solverChoice.coupling_type == CouplingType::TwoWay && finest_level > 0)
40  {
41  if (level < finest_level) {
42  fr_as_crse = getAdvFluxReg(level+1);
43  }
44  if (level > 0) {
45  fr_as_fine = getAdvFluxReg(level);
46  }
47  }
48 
49  MultiFab* new_detJ =
50  (solverChoice.terrain_type == TerrainType::MovingFittedMesh) ? (detJ_cc_new[level].get()) : (detJ_cc[level].get());
51 
52  erf_slow_rhs_post(level, finest_level, nrk, slow_dt, micro->Get_Qstate_Moist_Size(),
53  S_rhs, S_old, S_new, S_data, S_prim, avg_xmom[level], avg_ymom[level], avg_zmom[level],
54  xvel_new, yvel_new, zvel_new, cc_src, SmnSmn, eddyDiffs,
55  Hfx1, Hfx2, Hfx3, Q1fx1, Q1fx2, Q1fx3, Q2fx3, Diss,
56  fine_geom, solverChoice, m_SurfaceLayer, domain_bcs_type_d, domain_bcs_type,
57  z_phys_nd[level], z_phys_cc[level], ax[level], ay[level], az[level],
58  detJ_cc[level], new_detJ, stretched_dz_d[level], mapfac[level], EBFactory(level),
59 #if defined(ERF_USE_NETCDF)
60  moist_set_rhs,
61  start_time+old_stage_time,
62  start_bdy_time, final_bdy_time, bdy_time_interval,
63  real_width,
64  bdy_data_xlo, bdy_data_xhi, bdy_data_ylo, bdy_data_yhi,
65 #endif
66 #ifdef ERF_USE_SHOC
67  shoc_interface[level],
68 #endif
69  fr_as_crse, fr_as_fine,
70  m_r2d);
71 
72  // Implicit diffusion of moisture
73  #include "ERF_ImplicitPost.H"
74 
75  // Apply state redistribution for cons states
76 
77  if (solverChoice.terrain_type == TerrainType::EB)
78  {
79  Vector<int> is_valid_slow_var; is_valid_slow_var.resize(RhoQ1_comp+1,0);
80  if (solverChoice.turbChoice[level].use_tke) {is_valid_slow_var[ RhoKE_comp] = 1;}
81  is_valid_slow_var[RhoScalar_comp] = 1;
82  if (solverChoice.moisture_type != MoistureType::None) {
83  is_valid_slow_var[RhoQ1_comp] = 1;
84  }
85  const int num_comp_total = S_rhs[IntVars::cons].nComp();
86  const int num_grow = S_rhs[IntVars::cons].nGrow();
87  const int nvars = S_data[IntVars::cons].nComp();
88 
89  MultiFab dUdt_tmp(ba, dm, num_comp_total, num_grow, MFInfo(), EBFactory(level));
90  dUdt_tmp.setVal(0, 0, num_comp_total, num_grow);
91 
92  int start_comp;
93  int num_comp;
94 
95  for (int ivar(RhoKE_comp); ivar<= RhoQ1_comp; ++ivar)
96  {
97  if (is_valid_slow_var[ivar])
98  {
99  start_comp = ivar;
100  num_comp = 1;
101  if (ivar == RhoQ1_comp) {
102  num_comp = nvars - RhoQ1_comp;
103  } else if (ivar == RhoScalar_comp) {
104  num_comp = NSCALARS;
105  }
106  MultiFab::Copy(dUdt_tmp, S_rhs[IntVars::cons], start_comp, start_comp, num_comp, 0);
107  }
108  }
109  dUdt_tmp.FillBoundary(fine_geom.periodicity());
110  dUdt_tmp.setDomainBndry(amrex::Real(1.234e10), 0, num_comp_total, fine_geom);
111 
112  const BCRec* bc_ptr_d = domain_bcs_type_d.data();
113 
114  // Update S_rhs by Redistribution.
115  // To-do: Currently, redistributing all the scalar variables.
116  // This needs to be redistributed only for num_comp variables starting from ivar, for efficiency.
117  redistribute_term ( num_comp_total, fine_geom, S_rhs[IntVars::cons], dUdt_tmp,
118  S_old[IntVars::cons], EBFactory(level), bc_ptr_d, slow_dt);
119 
120  // Update state using the updated S_rhs. (NOTE: redistribute_term returns RHS not state variables.)
121  for ( MFIter mfi(S_new[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
122  {
123  Box tbx = mfi.tilebox();
124  const Array4<Real>& snew = S_new[IntVars::cons].array(mfi);
125  const Array4<Real>& sold = S_old[IntVars::cons].array(mfi);
126  const Array4<Real>& srhs = S_rhs[IntVars::cons].array(mfi);
127  Array4<const Real> detJ_arr = EBFactory(level).getVolFrac().const_array(mfi);
128 
129  for (int ivar(RhoKE_comp); ivar<= RhoQ1_comp; ++ivar)
130  {
131  if (is_valid_slow_var[ivar])
132  {
133  start_comp = ivar;
134  num_comp = 1;
135  if (ivar == RhoQ1_comp) {
136  num_comp = nvars - RhoQ1_comp;
137  } else if (ivar == RhoScalar_comp) {
138  num_comp = NSCALARS;
139  }
140  ParallelFor(tbx, num_comp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn)
141  {
142  if (detJ_arr(i,j,k) > zero) {
143  const int n = start_comp + nn;
144  snew(i,j,k,n) = sold(i,j,k,n) + slow_dt * srhs(i,j,k,n);
145  }
146  });
147  }
148  }
149  }
150  } // EB
151 
152  // Apply boundary conditions on all the state variables that have been updated
153  // in both the fast and slow integrators
154  apply_bcs(S_new, new_stage_time, S_new[IntVars::cons].nGrow(), S_new[IntVars::xmom].nGrow(),
155  fast_only=false, vel_and_mom_synced=false);
156 
157  if (solverChoice.moisture_tight_coupling) {
158  // TODO: need iteration var for lagrangian microphysics
159  // call signature in ERF::Advance() is
160  //advance_microphysics(lev, S_new, dt_lev, iteration, time);
161  advance_microphysics(level, S_new[0], slow_dt, 123456789, old_step_time);
162  }
163  }; // end slow_rhs_fun_post
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
@ nvars
Definition: ERF_DataStruct.H:98
void redistribute_term(int ncomp, const Geometry &geom, MultiFab &result, MultiFab &result_tmp, MultiFab const &state, EBFArrayBoxFactory const &ebfact, BCRec const *bc, Real const local_dt)
Definition: ERF_EBRedistribute.cpp:13
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
amrex::GpuArray< Real, AMREX_SPACEDIM > dxInv
Definition: ERF_InitCustomPertVels_ParticleTests.H:17
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);} } })
pp get("wavelength", wavelength)
amrex::Real Real
Definition: ERF_ShocInterface.H:19
void erf_slow_rhs_post(int level, int finest_level, int nrk, Real dt, int n_qstate, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_old, Vector< MultiFab > &S_new, Vector< MultiFab > &S_data, const MultiFab &S_prim, MultiFab &avg_xmom, MultiFab &avg_ymom, MultiFab &avg_zmom, const MultiFab &xvel, const MultiFab &yvel, const MultiFab &, const MultiFab &source, const MultiFab *SmnSmn, const MultiFab *eddyDiffs, MultiFab *Hfx1, MultiFab *Hfx2, MultiFab *Hfx3, MultiFab *Q1fx1, MultiFab *Q1fx2, MultiFab *Q1fx3, MultiFab *Q2fx3, MultiFab *Diss, const Geometry geom, const SolverChoice &solverChoice, std::unique_ptr< SurfaceLayer > &SurfLayer, const Gpu::DeviceVector< BCRec > &domain_bcs_type_d, const Vector< BCRec > &domain_bcs_type_h, std::unique_ptr< MultiFab > &z_phys_nd, std::unique_ptr< MultiFab > &z_phys_cc, std::unique_ptr< MultiFab > &ax, std::unique_ptr< MultiFab > &ay, std::unique_ptr< MultiFab > &az, std::unique_ptr< MultiFab > &detJ, MultiFab *detJ_new, Gpu::DeviceVector< Real > &stretched_dz_d, Vector< std::unique_ptr< MultiFab >> &mapfac, amrex::EBFArrayBoxFactory const &ebfact, YAFluxRegister *fr_as_crse, YAFluxRegister *fr_as_fine, std::unique_ptr< ReadBndryPlanes > &m_r2d)
Definition: ERF_SlowRhsPost.cpp:47
auto slow_rhs_fun_post
Definition: ERF_TI_slow_rhs_post.H:3
auto apply_bcs
Definition: ERF_TI_utils.H:73
@ cons
Definition: ERF_IndexDefines.H:176
@ xmom
Definition: ERF_IndexDefines.H:177