ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_TI_utils.H
Go to the documentation of this file.
1 /**
2  * Define the Exner pressure at this RK stage
3  */
4  auto make_pi_stage = [&](const MultiFab& cons_state)
5  {
6  BL_PROFILE("make_pi_stage()");
7 
8 #ifdef _OPENMP
9 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
10 #endif
11  for (MFIter mfi(cons_state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
12  {
13  //
14  // We only use one ghost cell of pi_stage so we only fill one here
15  //
16  const Box& gbx1 = mfi.growntilebox(1);
17 
18  const Array4<const Real>& cons_arr = cons_state.const_array(mfi);
19  const Array4< Real>& pi_stage_arr = pi_stage.array(mfi);
20  const Real rdOcp = solverChoice.rdOcp;
21  const bool l_use_moisture = (solverChoice.moisture_type != MoistureType::None);
22 
23  amrex::ParallelFor(gbx1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
24  {
25  auto qv_for_p = (l_use_moisture) ? cons_arr(i,j,k,RhoQ1_comp)/cons_arr(i,j,k,Rho_comp) : zero;
26  pi_stage_arr(i,j,k) = getExnergivenRTh(cons_arr(i,j,k,RhoTheta_comp), rdOcp, qv_for_p);
27  });
28  } // mfi
29  };
30 
31 /**
32  * This routine is called after the scalars are updated for each RK stage
33  */
34  auto apply_bcs = [&](Vector<MultiFab>& S_data,
35  const Real time_for_fp, int ng_cons, int ng_vel,
36  bool fast_only, bool vel_and_mom_synced)
37  {
38  BL_PROFILE("apply_bcs()");
39 
40  int scomp_cons;
41  int ncomp_cons;
42  bool cons_only;
43 
44  int ng_cons_to_use;
45 
46  // **********************************************************************************
47  // Because momentum is updated in the time-stepping routine, but boundary conditions
48  // are imposed on velocities, we must first update the velocity from the momentum
49  // before applying bcs.
50  // **********************************************************************************
51  if (!vel_and_mom_synced) {
52 
53  // **********************************************************************************
54  // Call FillPatch routines for the density only because we need it to convert between
55  // momentum and velocity
56  // This fills ghost cells/faces from
57  // 1) coarser level if lev > 0
58  // 2) physical boundaries
59  // 3) other grids at the same level
60  // **********************************************************************************
61 
62  // We must have at least one ghost cell of density to convert from momentum to velocity
63  // on the valid region
64  AMREX_ALWAYS_ASSERT (ng_cons >= 1);
65 
66  // We must have at least one extra ghost cell of density to convert from velocity to momentum
67  // on the valid region
68  ng_cons_to_use = std::max(ng_cons, ng_vel+1);
69 
70  scomp_cons = 0;
71  ncomp_cons = 1;
72  cons_only = true;
73  // **********************************************************************************
74  // Fill ghost cells of density only
75  // **********************************************************************************
76  FillIntermediatePatch(level, time_for_fp,
77  {&S_data[IntVars::cons], &xvel_new, &yvel_new, &zvel_new},
78  {&S_data[IntVars::cons], &S_data[IntVars::xmom],
79  &S_data[IntVars::ymom], &S_data[IntVars::zmom]},
80  ng_cons_to_use, 0, cons_only, scomp_cons, ncomp_cons);
81  }
82 
83  // ***************************************************************************************
84  // Call FillPatch routines for all data except rho which was filled above
85  // This fills ghost cells/faces from
86  // 1) coarser level if lev > 0
87  // 2) physical boundaries
88  // 3) other grids at the same level
89  // ***************************************************************************************
90  if (vel_and_mom_synced) {
91  if (fast_only) {
92  scomp_cons = 0;
93  ncomp_cons = 2; // rho and (rho theta) only
94  } else {
95  scomp_cons = 0;
96  ncomp_cons = S_data[IntVars::cons].nComp();
97  }
98  // We must have at least one extra ghost cell of density to convert from velocity to momentum
99  // on the valid region
100  ng_cons_to_use = std::max(ng_cons, ng_vel+1);
101 
102  } else {
103  if (fast_only) {
104  scomp_cons = 1;
105  ncomp_cons = 1; // (rho theta) only since we filled rho above
106  } else {
107  scomp_cons = 1;
108  ncomp_cons = S_data[IntVars::cons].nComp()-1; // since we filled rho above
109  }
110  ng_cons_to_use = ng_cons;
111  }
112 
113  // **********************************************************************************
114  // NOTE: FillIntermediatePatch takes momenta at the new time, and returns
115  // BOTH updated velocities and momenta
116  // **********************************************************************************
117  cons_only = false;
118  FillIntermediatePatch(level, time_for_fp,
119  {&S_data[IntVars::cons], &xvel_new, &yvel_new, &zvel_new},
120  {&S_data[IntVars::cons], &S_data[IntVars::xmom],
121  &S_data[IntVars::ymom], &S_data[IntVars::zmom]},
122  ng_cons_to_use, ng_vel, cons_only, scomp_cons, ncomp_cons);
123  };
124 
125  auto update_terrain_stage = [&](int lev, Real old_step_time, Real old_stage_time, Real new_stage_time, Real slow_dt)
126  {
127  // Note that the "old" and "new" metric terms correspond to
128  // t^n and the RK stage (either t^*, t^** or t^{n+1} that this source
129  // will be used to advance to
130 
131  // The "src" metric terms correspond to the time at which we are evaluating the source here,
132  // aka old_stage_time
133 
134  if (verbose) Print() << "Re-making old geometry at old time : " << old_step_time << std::endl;
135 
136  //
137  // Make a temporary MF to fill the terrain
138  //
139  Box terrain_bx(surroundingNodes(Geom(0).Domain())); terrain_bx.grow(z_phys_nd[lev]->nGrow());
140  FArrayBox terrain_fab_old(makeSlab(terrain_bx,2,0),1);
141  FArrayBox terrain_fab_src(makeSlab(terrain_bx,2,0),1);
142  FArrayBox terrain_fab_new(makeSlab(terrain_bx,2,0),1);
143 
144  prob->init_terrain_surface(fine_geom,terrain_fab_old,old_step_time);
145  prob->init_terrain_surface(fine_geom,terrain_fab_src,old_stage_time);
146  prob->init_terrain_surface(fine_geom,terrain_fab_new,new_stage_time);
147 
148  // Copy on intersection
149  for (MFIter mfi(*z_phys_nd[lev],false); mfi.isValid(); ++mfi)
150  {
151  Box isect = terrain_fab_old.box() & (*z_phys_nd[lev])[mfi].box();
152  ( *z_phys_nd[lev])[mfi].template copy<RunOn::Device>(terrain_fab_old,isect,0,isect,0,1);
153  (*z_phys_nd_src[lev])[mfi].template copy<RunOn::Device>(terrain_fab_src,isect,0,isect,0,1);
154  (*z_phys_nd_new[lev])[mfi].template copy<RunOn::Device>(terrain_fab_new,isect,0,isect,0,1);
155  }
156 
157  make_terrain_fitted_coords(lev,fine_geom,*z_phys_nd[lev], zlevels_stag[lev], phys_bc_type);
158  make_J (fine_geom, *z_phys_nd[lev], *detJ_cc[lev]);
159  make_areas (fine_geom, *z_phys_nd[lev], *ax[lev], *ay[lev], *az[lev]);
160 
161  if (verbose) Print() << "Making src geometry at old_stage_time: "
162  << std::setprecision(timeprecision) << old_stage_time << std::endl;
163  make_terrain_fitted_coords(lev,fine_geom,*z_phys_nd_src[lev], zlevels_stag[lev], phys_bc_type);
164  make_J (fine_geom, *z_phys_nd_src[lev], *detJ_cc_src[lev]);
165  make_areas (fine_geom, *z_phys_nd_src[lev], *ax_src[lev], *ay_src[lev], *az_src[lev]);
166  make_zcc (fine_geom, *z_phys_nd_src[lev], *z_phys_cc_src[lev]);
167 
168  if (verbose) Print() << "Making new geometry at new_stage_time: "
169  << std::setprecision(timeprecision) << new_stage_time << std::endl;
170  make_terrain_fitted_coords(lev,fine_geom,*z_phys_nd_new[lev], zlevels_stag[lev], phys_bc_type);
171  make_J (fine_geom, *z_phys_nd_new[lev], *detJ_cc_new[lev]);
172 
173  Real inv_dt = one/slow_dt;
174 
175 #ifdef _OPENMP
176 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
177 #endif
178  for (MFIter mfi(*z_t_rk[lev],TilingIfNotGPU()); mfi.isValid(); ++mfi)
179  {
180  Box gbx = mfi.growntilebox(IntVect(1,1,0));
181 
182  const Array4<Real >& z_t_arr = z_t_rk[lev]->array(mfi);
183  const Array4<Real const>& z_nd_new_arr = z_phys_nd_new[lev]->const_array(mfi);
184  const Array4<Real const>& z_nd_old_arr = z_phys_nd[lev]->const_array(mfi);
185 
186  // Loop over horizontal plane
187  amrex::ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
188  {
189  // Evaluate between RK stages assuming the geometry is linear between old and new time
190  z_t_arr(i,j,k) = fourth * inv_dt * (z_nd_new_arr(i+1,j+1,k) - z_nd_old_arr(i+1,j+1,k)
191  +z_nd_new_arr(i ,j+1,k) - z_nd_old_arr( i,j+1,k)
192  +z_nd_new_arr(i+1,j ,k) - z_nd_old_arr(i+1,j ,k)
193  +z_nd_new_arr(i ,j ,k) - z_nd_old_arr(i ,j ,k));
194  });
195 
196  } // mfi
197  };
198 
199  auto update_terrain_substep = [&](int lev, Real old_substep_time, Real new_substep_time, Real dtau,
200  Vector<MultiFab>& S_data, std::unique_ptr<MultiFab>& z_t_pert)
201  {
202  // Make "old" fast geom -- store in z_phys_nd for convenience
203  Box terrain_bx(surroundingNodes(Geom(0).Domain())); terrain_bx.grow(z_phys_nd[lev]->nGrow());
204 
205  if (verbose) Print() << "Making geometry at start of substep time: "
206  << std::setprecision(timeprecision) << old_substep_time << std::endl;
207  FArrayBox terrain_fab_old(makeSlab(terrain_bx,2,0),1);
208  prob->init_terrain_surface(fine_geom,terrain_fab_old,old_substep_time);
209 
210  // Make "new" fast geom
211  if (verbose) Print() << "Making geometry for end of substep time :"
212  << std::setprecision(timeprecision) << new_substep_time << std::endl;
213  FArrayBox terrain_fab_new(makeSlab(terrain_bx,2,0),1);
214  prob->init_terrain_surface(fine_geom,terrain_fab_new,new_substep_time);
215 
216  // Copy on intersection
217  for (MFIter mfi(*z_phys_nd[lev],false); mfi.isValid(); ++mfi)
218  {
219  Box isect = terrain_fab_old.box() & (*z_phys_nd[lev])[mfi].box();
220  ( *z_phys_nd[lev])[mfi].template copy<RunOn::Device>(terrain_fab_old,isect,0,isect,0,1);
221  (*z_phys_nd_new[lev])[mfi].template copy<RunOn::Device>(terrain_fab_new,isect,0,isect,0,1);
222  }
223 
224  make_terrain_fitted_coords(lev,fine_geom,*z_phys_nd[lev], zlevels_stag[lev], phys_bc_type);
225  make_J (fine_geom,*z_phys_nd[lev], *detJ_cc[lev]);
226 
227  make_terrain_fitted_coords(lev,fine_geom,*z_phys_nd_new[lev], zlevels_stag[lev], phys_bc_type);
228  make_J (fine_geom,*z_phys_nd_new[lev], *detJ_cc_new[lev]);
229 
230  Real inv_dt = one/dtau;
231 
232  z_t_pert = std::make_unique<MultiFab>(S_data[IntVars::zmom].boxArray(), S_data[IntVars::zmom].DistributionMap(), 1, 1);
233 
234  for (MFIter mfi(*z_t_rk[lev],TilingIfNotGPU()); mfi.isValid(); ++mfi)
235  {
236  Box gbx = mfi.growntilebox(IntVect(1,1,0));
237 
238  const Array4<Real >& z_t_arr = z_t_rk[lev]->array(mfi);
239  const Array4<Real >& zp_t_arr = z_t_pert->array(mfi);
240 
241  const Array4<Real const>& z_nd_new_arr = z_phys_nd_new[lev]->const_array(mfi);
242  const Array4<Real const>& z_nd_old_arr = z_phys_nd[lev]->const_array(mfi);
243 
244  // Loop over horizontal plane
245  amrex::ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
246  {
247  // Evaluate between RK stages assuming the geometry is linear between old and new time
248  zp_t_arr(i,j,k) = fourth * inv_dt * (z_nd_new_arr(i+1,j+1,k) - z_nd_old_arr(i+1,j+1,k)
249  +z_nd_new_arr(i ,j+1,k) - z_nd_old_arr( i,j+1,k)
250  +z_nd_new_arr(i+1,j ,k) - z_nd_old_arr(i+1,j ,k)
251  +z_nd_new_arr(i ,j ,k) - z_nd_old_arr(i ,j ,k));
252  // Convert to perturbation: z"_t(t) = z_t(t) - z_t^{RK}
253  zp_t_arr(i,j,k) -= z_t_arr(i,j,k);
254  });
255  } // mfi
256  };
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenRTh(const amrex::Real rhotheta, const amrex::Real rdOcp, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:156
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
const Real rdOcp
Definition: ERF_InitCustomPert_Bomex.H:16
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
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
auto make_pi_stage
Definition: ERF_TI_utils.H:4
auto apply_bcs
Definition: ERF_TI_utils.H:34
auto update_terrain_substep
Definition: ERF_TI_utils.H:199
auto update_terrain_stage
Definition: ERF_TI_utils.H:125
void make_areas(const Geometry &geom, MultiFab &z_phys_nd, MultiFab &ax, MultiFab &ay, MultiFab &az)
Definition: ERF_TerrainMetrics.cpp:561
void make_terrain_fitted_coords(int lev, const Geometry &geom, MultiFab &z_phys_nd, Vector< Real > const &z_levels_h, GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
Definition: ERF_TerrainMetrics.cpp:47
void make_J(const Geometry &geom, MultiFab &z_phys_nd, MultiFab &detJ_cc)
Definition: ERF_TerrainMetrics.cpp:523
void make_zcc(const Geometry &geom, MultiFab &z_phys_nd, MultiFab &z_phys_cc)
Definition: ERF_TerrainMetrics.cpp:628
@ ymom
Definition: ERF_IndexDefines.H:194
@ cons
Definition: ERF_IndexDefines.H:192
@ zmom
Definition: ERF_IndexDefines.H:195
@ xmom
Definition: ERF_IndexDefines.H:193