52 template <
typename MF>
53 std::enable_if_t<IsFabArray<MF>::value>
operator() (
int, MF&)
const {}
60void rk_update (MF& Unew, MF
const& Uold, MF
const& dUdt,
Real dt)
62 auto const& snew = Unew.arrays();
63 auto const& sold = Uold.const_arrays();
64 auto const& sdot = dUdt.const_arrays();
66 (
int bi,
int i,
int j,
int k,
int n)
noexcept
68 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*sdot[bi](i,j,k,n);
70 Gpu::streamSynchronize();
75void rk_update (MF& Unew, MF
const& Uold, MF
const& dUdt1, MF
const& dUdt2, Real dt)
77 auto const& snew = Unew.arrays();
78 auto const& sold = Uold.const_arrays();
79 auto const& sdot1 = dUdt1.const_arrays();
80 auto const& sdot2 = dUdt2.const_arrays();
82 (
int bi,
int i,
int j,
int k,
int n)
noexcept
84 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*(sdot1[bi](i,j,k,n) +
87 Gpu::streamSynchronize();
92void rk2_update_2 (MF& Unew, MF
const& Uold, MF
const& dUdt, Real dt)
94 auto const& snew = Unew.arrays();
95 auto const& sold = Uold.const_arrays();
96 auto const& sdot = dUdt.const_arrays();
98 (
int bi,
int i,
int j,
int k,
int n)
noexcept
100 snew[bi](i,j,k,n) =
Real(0.5)*(snew[bi](i,j,k,n) +
102 sdot[bi](i,j,k,n) * dt);
104 Gpu::streamSynchronize();
108template <
typename MF>
109void rk3_update_3 (MF& Unew, MF
const& Uold, Array<MF,3>
const& rkk, Real dt6)
111 auto const& snew = Unew.arrays();
112 auto const& sold = Uold.const_arrays();
113 auto const& k1 = rkk[0].const_arrays();
114 auto const& k2 = rkk[1].const_arrays();
115 auto const& k3 = rkk[2].const_arrays();
117 (
int bi,
int i,
int j,
int k,
int n)
noexcept
119 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
120 + dt6 * (k1[bi](i,j,k,n) + k2[bi](i,j,k,n)
121 +
Real(4.) * k3[bi](i,j,k,n));
123 Gpu::streamSynchronize();
127template <
typename MF>
128void rk4_update_4 (MF& Unew, MF
const& Uold, Array<MF,4>
const& rkk, Real dt6)
130 auto const& snew = Unew.arrays();
131 auto const& sold = Uold.const_arrays();
132 auto const& k1 = rkk[0].const_arrays();
133 auto const& k2 = rkk[1].const_arrays();
134 auto const& k3 = rkk[2].const_arrays();
135 auto const& k4 = rkk[3].const_arrays();
137 (
int bi,
int i,
int j,
int k,
int n)
noexcept
139 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
140 + dt6 * ( k1[bi](i,j,k,n) + k4[bi](i,j,k,n)
141 +
Real(2.)*(k2[bi](i,j,k,n) + k3[bi](i,j,k,n)));
143 Gpu::streamSynchronize();
159template <
typename MF,
typename F,
typename FB,
typename P = PostStageNoOp>
160void RK2 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
165 MF dUdt(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
166 MFInfo(), Unew.Factory());
169 fillbndry(1, Uold, time);
170 frhs(1, dUdt, Uold, time,
Real(0.5)*dt);
172 detail::rk_update(Unew, Uold, dUdt, dt);
176 fillbndry(2, Unew, time+dt);
177 frhs(2, dUdt, Unew, time,
Real(0.5)*dt);
180 detail::rk2_update_2(Unew, Uold, dUdt, dt);
196template <
typename MF,
typename F,
typename FB,
typename R,
197 typename P = PostStageNoOp>
198void RK3 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
199 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
204 for (
auto& mf : rkk) {
205 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
206 MFInfo(), Unew.Factory());
210 fillbndry(1, Uold, time);
211 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
213 detail::rk_update(Unew, Uold, rkk[0], dt);
217 fillbndry(2, Unew, time+dt);
218 frhs(2, rkk[1], Unew, time+dt, dt/
Real(6.));
220 detail::rk_update(Unew, Uold, rkk[0], rkk[1],
Real(0.25)*dt);
225 fillbndry(3, Unew, t_half);
226 frhs(3, rkk[2], Unew, t_half, dt*
Real(2./3.));
228 detail::rk3_update_3(Unew, Uold, rkk,
Real(1./6.)*dt);
231 store_crse_data(rkk);
246template <
typename MF,
typename F,
typename FB,
typename R,
247 typename P = PostStageNoOp>
248void RK4 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
249 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
254 for (
auto& mf : rkk) {
255 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
256 MFInfo(), Unew.Factory());
260 fillbndry(1, Uold, time);
261 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
263 detail::rk_update(Unew, Uold, rkk[0],
Real(0.5)*dt);
268 fillbndry(2, Unew, t_half);
269 frhs(2, rkk[1], Unew, t_half, dt/
Real(3.));
271 detail::rk_update(Unew, Uold, rkk[1],
Real(0.5)*dt);
275 fillbndry(3, Unew, t_half);
276 frhs(3, rkk[2], Unew, t_half, dt/
Real(3.));
278 detail::rk_update(Unew, Uold, rkk[2], dt);
282 fillbndry(4, Unew, time+dt);
283 frhs(4, rkk[3], Unew, time+dt, dt/
Real(6.));
285 detail::rk4_update_4(Unew, Uold, rkk,
Real(1./6.)*dt);
288 store_crse_data(rkk);
void RK2(MF &Uold, MF &Unew, Real time, Real dt, F const &frhs, FB const &fillbndry, P const &post_stage=PostStageNoOp())
Time stepping with RK2.
Definition AMReX_RungeKutta.H:160
void RK3(MF &Uold, MF &Unew, Real time, Real dt, F const &frhs, FB const &fillbndry, R const &store_crse_data, P const &post_stage=PostStageNoOp())
Time stepping with RK3.
Definition AMReX_RungeKutta.H:198
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:193
Definition AMReX_RungeKutta.H:51
std::enable_if_t< IsFabArray< MF >::value > operator()(int, MF &) const
Definition AMReX_RungeKutta.H:53