5 Vector<MultiFab>& S_old,
6 Vector<MultiFab>& F_slow,
7 const Real time_for_fp,
const Real slow_dt,
10 BL_PROFILE(
"no_substep_fun");
11 amrex::ignore_unused(nrk);
14 const auto& dxInv = fine_geom.InvCellSizeArray();
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};
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;
26 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
29 for ( MFIter mfi(S_sum[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
31 const Box bx = mfi.tilebox();
32 Box tbx = mfi.nodaltilebox(0);
33 Box tby = mfi.nodaltilebox(1);
34 Box tbz = mfi.nodaltilebox(2);
36 Vector<Array4<Real> > ssum_h(n_data);
37 Vector<Array4<Real> > sold_h(n_data);
38 Vector<Array4<Real> > fslow_h(n_data);
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);
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);
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());
54 Array4<Real>* sold = sold_d.dataPtr();
55 Array4<Real>* ssum = ssum_d.dataPtr();
56 Array4<Real>* fslow = fslow_d.dataPtr();
59 if ( solverChoice.terrain_type == TerrainType::MovingFittedMesh )
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);
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);
67 const Array4<Real >& z_t_arr = z_t_rk[level]->array(mfi);
71 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int nn) {
80 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
86 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
93 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
102 Real dJ_old_kface = 0.5 * (dJ_old(i,j,k) + dJ_old(i,j,k-1));
103 Real dJ_new_kface = 0.5 * (dJ_new(i,j,k) + dJ_new(i,j,k-1));
110 const Array4<const Real>& dJ_old = detJ_cc[level]->const_array(mfi);
112 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int nn) {
114 if (dJ_old(i,j,k) > 0.0) {
121 if (solverChoice.terrain_type == TerrainType::EB)
123 ParallelFor(tbx, tby, tbz,
124 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
127 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
130 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
135 ParallelFor(tbx, tby, tbz,
136 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
140 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
144 [=] AMREX_GPU_DEVICE (
int i,
int j,
int k) noexcept {
152 if (solverChoice.terrain_type == TerrainType::EB)
159 dUdt_tmp.FillBoundary(fine_geom.periodicity());
160 dUdt_tmp.setDomainBndry(1.234e10, 0, ncomp_fast[
IntVars::cons], fine_geom);
162 BCRec
const* bc_ptr_d = domain_bcs_type_d.data();
172 for ( MFIter mfi(S_sum[
IntVars::cons],TilingIfNotGPU()); mfi.isValid(); ++mfi)
174 const Box bx = mfi.tilebox();
175 Vector<Array4<Real> > ssum_h(n_data);
176 Vector<Array4<Real> > sold_h(n_data);
177 Vector<Array4<Real> > fslow_h(n_data);
178 const Array4<const Real>& dJ_old = detJ_cc[level]->const_array(mfi);
180 for (
int i = 0; i < n_data; ++i) {
181 ssum_h[i] = S_sum[i].array(mfi);
182 sold_h[i] = S_old[i].array(mfi);
183 fslow_h[i] = F_slow[i].array(mfi);
186 Gpu::AsyncVector<Array4<Real> > sold_d(n_data);
187 Gpu::AsyncVector<Array4<Real> > ssum_d(n_data);
188 Gpu::AsyncVector<Array4<Real> > fslow_d(n_data);
190 Gpu::copy(Gpu::hostToDevice, sold_h.begin(), sold_h.end(), sold_d.begin());
191 Gpu::copy(Gpu::hostToDevice, ssum_h.begin(), ssum_h.end(), ssum_d.begin());
192 Gpu::copy(Gpu::hostToDevice, fslow_h.begin(), fslow_h.end(), fslow_d.begin());
194 Array4<Real>* sold = sold_d.dataPtr();
195 Array4<Real>* ssum = ssum_d.dataPtr();
196 Array4<Real>* fslow = fslow_d.dataPtr();
198 ParallelFor(bx, ncomp_fast[
IntVars::cons], [=] AMREX_GPU_DEVICE (
int i,
int j,
int k,
int nn)
201 if (dJ_old(i,j,k) > 0.0) {
215 fast_only=
true, vel_and_mom_synced=
false);
217 if (solverChoice.anelastic[level]) {
218 bool have_tb = (thin_xforce[0] || thin_yforce[0] || thin_zforce[0]);
220 project_velocities_tb(level, slow_dt, S_sum, pp_inc[level]);
222 project_velocities(level, slow_dt, S_sum, pp_inc[level]);
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:11
#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:50
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int i, int j, int k, amrex::Real omega, amrex::Real u, amrex::Real v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:414
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:94
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:137
@ NumTypes
Definition: ERF_IndexDefines.H:154
@ ymom
Definition: ERF_IndexDefines.H:152
@ cons
Definition: ERF_IndexDefines.H:150
@ zmom
Definition: ERF_IndexDefines.H:153
@ xmom
Definition: ERF_IndexDefines.H:151