ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_TI_no_substep_fun.H
Go to the documentation of this file.
1 /**
2  * Wrapper for advancing the solution with the slow RHS in the absence of acoustic substepping
3  */
4  auto no_substep_fun = [&](Vector<MultiFab>& S_sum,
5  Vector<MultiFab>& S_old,
6  Vector<MultiFab>& F_slow,
7  const Real time_for_fp, const Real slow_dt,
8  const int nrk)
9  {
10  BL_PROFILE("no_substep_fun");
11  amrex::ignore_unused(nrk);
12  int n_data = IntVars::NumTypes;
13 
14  const auto& dxInv = fine_geom.InvCellSizeArray();
15 
16  const amrex::GpuArray<int, IntVars::NumTypes> scomp_fast = {0,0,0,0};
17  const amrex::GpuArray<int, IntVars::NumTypes> ncomp_fast = {2,1,1,1};
18 
19  if (verbose) amrex::Print() << " No-substepping time integration at level " << level
20  << std::setprecision(timeprecision)
21  << " to " << time_for_fp
22  << " with dt = " << slow_dt << std::endl;
23 
24  // Update S_sum = S_stage only for the fast variables
25 #ifdef _OPENMP
26 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
27 #endif
28  {
29  for ( MFIter mfi(S_sum[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
30  {
31  const Box bx = mfi.tilebox();
32  Box tbx = mfi.nodaltilebox(0);
33  Box tby = mfi.nodaltilebox(1);
34  Box tbz = mfi.nodaltilebox(2);
35 
36  Vector<Array4<Real> > ssum_h(n_data);
37  Vector<Array4<Real> > sold_h(n_data);
38  Vector<Array4<Real> > fslow_h(n_data);
39 
40  for (int i = 0; i < n_data; ++i) {
41  ssum_h[i] = S_sum[i].array(mfi);
42  sold_h[i] = S_old[i].array(mfi);
43  fslow_h[i] = F_slow[i].array(mfi);
44  }
45 
46  Gpu::AsyncVector<Array4<Real> > sold_d(n_data);
47  Gpu::AsyncVector<Array4<Real> > ssum_d(n_data);
48  Gpu::AsyncVector<Array4<Real> > fslow_d(n_data);
49 
50  Gpu::copy(Gpu::hostToDevice, sold_h.begin(), sold_h.end(), sold_d.begin());
51  Gpu::copy(Gpu::hostToDevice, ssum_h.begin(), ssum_h.end(), ssum_d.begin());
52  Gpu::copy(Gpu::hostToDevice, fslow_h.begin(), fslow_h.end(), fslow_d.begin());
53 
54  Array4<Real>* sold = sold_d.dataPtr();
55  Array4<Real>* ssum = ssum_d.dataPtr();
56  Array4<Real>* fslow = fslow_d.dataPtr();
57 
58  // Moving terrain
59  if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
60  {
61  const Array4<const Real>& dJ_old = detJ_cc[level]->const_array(mfi);
62  const Array4<const Real>& dJ_new = detJ_cc_new[level]->const_array(mfi);
63 
64  const Array4<const Real>& z_nd_old = z_phys_nd[level]->const_array(mfi);
65  const Array4<const Real>& z_nd_new = z_phys_nd_new[level]->const_array(mfi);
66 
67  const Array4<const Real>& mf_u = mapfac[level][MapFacType::u_x]->const_array(mfi);
68  const Array4<const Real>& mf_v = mapfac[level][MapFacType::v_y]->const_array(mfi);
69 
70  const Array4<Real >& z_t_arr = z_t_rk[level]->array(mfi);
71 
72  // We have already scaled the slow source term to have the extra factor of dJ
73  ParallelFor(bx, ncomp_fast[IntVars::cons],
74  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) {
75  const int n = scomp_fast[IntVars::cons] + nn;
76  ssum[IntVars::cons](i,j,k,n) = dJ_old(i,j,k) * sold[IntVars::cons](i,j,k,n)
77  + slow_dt * fslow[IntVars::cons](i,j,k,n);
78  ssum[IntVars::cons](i,j,k,n) /= dJ_new(i,j,k);
79  });
80 
81  // We have already scaled the slow source term to have the extra factor of dJ
82  ParallelFor(tbx, tby,
83  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
84  Real h_zeta_old = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd_old);
85  Real h_zeta_new = Compute_h_zeta_AtIface(i, j, k, dxInv, z_nd_new);
86  ssum[IntVars::xmom](i,j,k) = ( h_zeta_old * sold[IntVars::xmom](i,j,k)
87  + slow_dt * fslow[IntVars::xmom](i,j,k) ) / h_zeta_new;
88  },
89  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
90  Real h_zeta_old = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd_old);
91  Real h_zeta_new = Compute_h_zeta_AtJface(i, j, k, dxInv, z_nd_new);
92  ssum[IntVars::ymom](i,j,k) = ( h_zeta_old * sold[IntVars::ymom](i,j,k)
93  + slow_dt * fslow[IntVars::ymom](i,j,k) ) / h_zeta_new;
94  });
95  ParallelFor(tbz,
96  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
97  if (k == 0) {
98  // Here we take advantage of the fact that moving terrain has a slip wall
99  // so we can just use the new value at (i,j,0).
100  Real rho_on_face = ssum[IntVars::cons](i,j,k,Rho_comp);
101  ssum[IntVars::zmom](i,j,k) = WFromOmega(i,j,k,rho_on_face*z_t_arr(i,j,k),
102  ssum[IntVars::xmom], ssum[IntVars::ymom],
103  mf_u, mf_v, z_nd_new,dxInv);
104  } else {
105  Real dJ_old_kface = 0.5 * (dJ_old(i,j,k) + dJ_old(i,j,k-1));
106  Real dJ_new_kface = 0.5 * (dJ_new(i,j,k) + dJ_new(i,j,k-1));
107  ssum[IntVars::zmom](i,j,k) = ( dJ_old_kface * sold[IntVars::zmom](i,j,k)
108  + slow_dt * fslow[IntVars::zmom](i,j,k) ) / dJ_new_kface;
109  }
110  });
111 
112  } else { // Fixed or no terrain
113  const Array4<const Real>& dJ_old = detJ_cc[level]->const_array(mfi);
114  ParallelFor(bx, ncomp_fast[IntVars::cons],
115  [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn) {
116  const int n = scomp_fast[IntVars::cons] + nn;
117  if (dJ_old(i,j,k) > 0.0) {
118  ssum[IntVars::cons](i,j,k,n) = sold[IntVars::cons](i,j,k,n) + slow_dt *
119  ( fslow[IntVars::cons](i,j,k,n) );
120  }
121  });
122 
123  // Commenting out the update is a HACK while developing the EB capability
124  if (solverChoice.terrain_type == TerrainType::EB)
125  {
126  const Array4<const Real>& vfrac_u = (get_eb(level).get_u_const_factory())->getVolFrac().const_array(mfi);
127  const Array4<const Real>& vfrac_v = (get_eb(level).get_v_const_factory())->getVolFrac().const_array(mfi);
128  const Array4<const Real>& vfrac_w = (get_eb(level).get_w_const_factory())->getVolFrac().const_array(mfi);
129 
130  ParallelFor(tbx, tby, tbz,
131  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
132  if (vfrac_u(i,j,k) > 0.0) {
133  // ssum[IntVars::xmom](i,j,k) = sold[IntVars::xmom](i,j,k);
134  ssum[IntVars::xmom](i,j,k) = sold[IntVars::xmom](i,j,k)
135  + slow_dt * fslow[IntVars::xmom](i,j,k);
136  } else {
137  ssum[IntVars::xmom](i,j,k) = sold[IntVars::xmom](i,j,k);
138  }
139  },
140  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
141  if (vfrac_v(i,j,k) > 0.0) {
142  // ssum[IntVars::ymom](i,j,k) = sold[IntVars::ymom](i,j,k);
143  ssum[IntVars::ymom](i,j,k) = sold[IntVars::ymom](i,j,k)
144  + slow_dt * fslow[IntVars::ymom](i,j,k);
145  } else {
146  ssum[IntVars::ymom](i,j,k) = sold[IntVars::ymom](i,j,k);
147  }
148  },
149  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
150  if (vfrac_w(i,j,k) > 0.0) {
151  // ssum[IntVars::zmom](i,j,k) = sold[IntVars::zmom](i,j,k);
152  ssum[IntVars::zmom](i,j,k) = sold[IntVars::zmom](i,j,k)
153  + slow_dt * fslow[IntVars::zmom](i,j,k);
154  } else {
155  ssum[IntVars::zmom](i,j,k) = sold[IntVars::zmom](i,j,k);
156  }
157  });
158 
159  } else {
160  ParallelFor(tbx, tby, tbz,
161  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
162  ssum[IntVars::xmom](i,j,k) = sold[IntVars::xmom](i,j,k)
163  + slow_dt * fslow[IntVars::xmom](i,j,k);
164  },
165  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
166  ssum[IntVars::ymom](i,j,k) = sold[IntVars::ymom](i,j,k)
167  + slow_dt * fslow[IntVars::ymom](i,j,k);
168  },
169  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
170  ssum[IntVars::zmom](i,j,k) = sold[IntVars::zmom](i,j,k)
171  + slow_dt * fslow[IntVars::zmom](i,j,k);
172  });
173  } // no EB
174  } // not moving terrain
175  } // mfi
176 
177  if (solverChoice.terrain_type == TerrainType::EB)
178  {
179  // Redistribute cons states (cell-centered)
180 
181  MultiFab dUdt_tmp(ba, dm, F_slow[IntVars::cons].nComp(), F_slow[IntVars::cons].nGrow(), MFInfo(), EBFactory(level));
182  dUdt_tmp.setVal(0.);
183  MultiFab::Copy(dUdt_tmp, F_slow[IntVars::cons], 0, 0, F_slow[IntVars::cons].nComp(), F_slow[IntVars::cons].nGrow());
184  dUdt_tmp.FillBoundary(fine_geom.periodicity());
185  dUdt_tmp.setDomainBndry(1.234e10, 0, ncomp_fast[IntVars::cons], fine_geom);
186 
187  BCRec const* bc_ptr_d = domain_bcs_type_d.data();
188  S_old[IntVars::cons].FillBoundary(fine_geom.periodicity());
189  S_old[IntVars::cons].setDomainBndry(1.234e10, 0, ncomp_fast[IntVars::cons], fine_geom);
190 
191  // Update F_slow by Redistribution.
192  redistribute_term ( ncomp_fast[IntVars::cons], fine_geom, F_slow[IntVars::cons], dUdt_tmp,
193  S_old[IntVars::cons], EBFactory(level), bc_ptr_d, slow_dt);
194 
195  // Redistribute momentum states (face-centered)
196 
197  Vector<MultiFab> dUdt_mom(AMREX_SPACEDIM);
198  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
199  dUdt_mom[i].define(F_slow[1+i].boxArray(), F_slow[1+i].DistributionMap(), F_slow[1+i].nComp(), F_slow[1+i].nGrow(), MFInfo());
200  dUdt_mom[i].setVal(0.);
201  MultiFab::Copy(dUdt_mom[i], F_slow[1+i], 0, 0, F_slow[1+i].nComp(), F_slow[1+i].nGrow());
202  dUdt_mom[i].FillBoundary(fine_geom.periodicity());
203  dUdt_mom[i].setDomainBndry(1.234e10, 0, ncomp_fast[1+i], fine_geom);
204  }
205 
206  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
207  S_old[1+i].FillBoundary(fine_geom.periodicity());
208  S_old[1+i].setDomainBndry(1.234e10, 0, ncomp_fast[1+i], fine_geom);
209  }
210 
211  // Update F_slow by Redistribution.
212  redistribute_term( ncomp_fast[IntVars::xmom], fine_geom, F_slow[IntVars::xmom], dUdt_mom[0],
213  S_old[IntVars::xmom], *(get_eb(level).get_u_const_factory()), bc_ptr_d, slow_dt);
214  redistribute_term( ncomp_fast[IntVars::ymom], fine_geom, F_slow[IntVars::ymom], dUdt_mom[1],
215  S_old[IntVars::ymom], *(get_eb(level).get_v_const_factory()), bc_ptr_d, slow_dt);
216  redistribute_term( ncomp_fast[IntVars::zmom], fine_geom, F_slow[IntVars::zmom], dUdt_mom[2],
217  S_old[IntVars::zmom], *(get_eb(level).get_w_const_factory()), bc_ptr_d, slow_dt);
218 
219  // Update state using the updated F_slow. (NOTE: redistribute_term returns RHS not state variables.)
220 
221  for ( MFIter mfi(S_sum[IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
222  {
223  const Box bx = mfi.tilebox();
224  Box tbx = mfi.nodaltilebox(0);
225  Box tby = mfi.nodaltilebox(1);
226  Box tbz = mfi.nodaltilebox(2);
227 
228  Vector<Array4<Real> > ssum_h(n_data);
229  Vector<Array4<Real> > sold_h(n_data);
230  Vector<Array4<Real> > fslow_h(n_data);
231 
232  for (int i = 0; i < n_data; ++i) {
233  ssum_h[i] = S_sum[i].array(mfi);
234  sold_h[i] = S_old[i].array(mfi);
235  fslow_h[i] = F_slow[i].array(mfi);
236  }
237 
238  Gpu::AsyncVector<Array4<Real> > sold_d(n_data);
239  Gpu::AsyncVector<Array4<Real> > ssum_d(n_data);
240  Gpu::AsyncVector<Array4<Real> > fslow_d(n_data);
241 
242  Gpu::copy(Gpu::hostToDevice, sold_h.begin(), sold_h.end(), sold_d.begin());
243  Gpu::copy(Gpu::hostToDevice, ssum_h.begin(), ssum_h.end(), ssum_d.begin());
244  Gpu::copy(Gpu::hostToDevice, fslow_h.begin(), fslow_h.end(), fslow_d.begin());
245 
246  Array4<Real>* sold = sold_d.dataPtr();
247  Array4<Real>* ssum = ssum_d.dataPtr();
248  Array4<Real>* fslow = fslow_d.dataPtr();
249 
250  // const Array4<const Real>& vfrac_c = detJ_cc[level]->const_array(mfi);
251  const Array4<const Real>& vfrac_c = (get_eb(level).get_const_factory())->getVolFrac().const_array(mfi);
252 
253  ParallelFor(bx, ncomp_fast[IntVars::cons], [=] AMREX_GPU_DEVICE (int i, int j, int k, int nn)
254  {
255  const int n = scomp_fast[IntVars::cons] + nn;
256  if (vfrac_c(i,j,k) > 0.0) {
257  ssum[IntVars::cons](i,j,k,n) = sold[IntVars::cons](i,j,k,n) + slow_dt *
258  ( fslow[IntVars::cons](i,j,k,n) );
259  }
260  });
261 
262  const Array4<const Real>& vfrac_u = (get_eb(level).get_u_const_factory())->getVolFrac().const_array(mfi);
263  const Array4<const Real>& vfrac_v = (get_eb(level).get_v_const_factory())->getVolFrac().const_array(mfi);
264  const Array4<const Real>& vfrac_w = (get_eb(level).get_w_const_factory())->getVolFrac().const_array(mfi);
265 
266  ParallelFor(tbx, tby, tbz,
267  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
268  if (vfrac_u(i,j,k) > 0.0) {
269  ssum[IntVars::xmom](i,j,k) = sold[IntVars::xmom](i,j,k)
270  + slow_dt * fslow[IntVars::xmom](i,j,k);
271  }
272  },
273  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
274  if (vfrac_v(i,j,k) > 0.0) {
275  ssum[IntVars::ymom](i,j,k) = sold[IntVars::ymom](i,j,k)
276  + slow_dt * fslow[IntVars::ymom](i,j,k);
277  }
278  },
279  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
280  if (vfrac_w(i,j,k) > 0.0) {
281  ssum[IntVars::zmom](i,j,k) = sold[IntVars::zmom](i,j,k)
282  + slow_dt * fslow[IntVars::zmom](i,j,k);
283  }
284  });
285  } // MFIter
286 
287  } // EB
288 
289  } // omp
290 
291  // Even if we update all the conserved variables we don't need
292  // to fillpatch the slow ones every acoustic substep
293  apply_bcs(S_sum, time_for_fp, S_sum[IntVars::cons].nGrow(), S_sum[IntVars::xmom].nGrow(),
294  fast_only=true, vel_and_mom_synced=false);
295 
296  if (solverChoice.anelastic[level]) {
297  bool have_tb = (thin_xforce[0] || thin_yforce[0] || thin_zforce[0]);
298  if (have_tb) {
299  project_velocity_tb(level, slow_dt, S_sum);
300  } else {
301  project_momenta(level, slow_dt, S_sum);
302  }
303  }
304  };
@ v_y
Definition: ERF_DataStruct.H:23
@ u_x
Definition: ERF_DataStruct.H:22
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 Rho_comp
Definition: ERF_IndexDefines.H:36
auto no_substep_fun
Definition: ERF_TI_no_substep_fun.H:4
auto apply_bcs
Definition: ERF_TI_utils.H:71
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:96
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:139
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int &i, int &j, int &k, amrex::Real omega, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:465
@ NumTypes
Definition: ERF_IndexDefines.H:162
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159