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 if (!Gpu::inNoSyncRegion()) {
71 Gpu::streamSynchronize();
77void rk_update (MF& Unew, MF
const& Uold, MF
const& dUdt1, MF
const& dUdt2, Real dt)
79 auto const& snew = Unew.arrays();
80 auto const& sold = Uold.const_arrays();
81 auto const& sdot1 = dUdt1.const_arrays();
82 auto const& sdot2 = dUdt2.const_arrays();
84 (
int bi,
int i,
int j,
int k,
int n)
noexcept
86 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*(sdot1[bi](i,j,k,n) +
89 if (!Gpu::inNoSyncRegion()) {
90 Gpu::streamSynchronize();
96void rk2_update_2 (MF& Unew, MF
const& Uold, MF
const& dUdt, Real dt)
98 auto const& snew = Unew.arrays();
99 auto const& sold = Uold.const_arrays();
100 auto const& sdot = dUdt.const_arrays();
102 (
int bi,
int i,
int j,
int k,
int n)
noexcept
104 snew[bi](i,j,k,n) =
Real(0.5)*(snew[bi](i,j,k,n) +
106 sdot[bi](i,j,k,n) * dt);
108 if (!Gpu::inNoSyncRegion()) {
109 Gpu::streamSynchronize();
114template <
typename MF>
115void rk3_update_3 (MF& Unew, MF
const& Uold, Array<MF,3>
const& rkk, Real dt6)
117 auto const& snew = Unew.arrays();
118 auto const& sold = Uold.const_arrays();
119 auto const& k1 = rkk[0].const_arrays();
120 auto const& k2 = rkk[1].const_arrays();
121 auto const& k3 = rkk[2].const_arrays();
123 (
int bi,
int i,
int j,
int k,
int n)
noexcept
125 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
126 + dt6 * (k1[bi](i,j,k,n) + k2[bi](i,j,k,n)
127 +
Real(4.) * k3[bi](i,j,k,n));
129 if (!Gpu::inNoSyncRegion()) {
130 Gpu::streamSynchronize();
135template <
typename MF>
136void rk4_update_4 (MF& Unew, MF
const& Uold, Array<MF,4>
const& rkk, Real dt6)
138 auto const& snew = Unew.arrays();
139 auto const& sold = Uold.const_arrays();
140 auto const& k1 = rkk[0].const_arrays();
141 auto const& k2 = rkk[1].const_arrays();
142 auto const& k3 = rkk[2].const_arrays();
143 auto const& k4 = rkk[3].const_arrays();
145 (
int bi,
int i,
int j,
int k,
int n)
noexcept
147 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
148 + dt6 * ( k1[bi](i,j,k,n) + k4[bi](i,j,k,n)
149 +
Real(2.)*(k2[bi](i,j,k,n) + k3[bi](i,j,k,n)));
151 if (!Gpu::inNoSyncRegion()) {
152 Gpu::streamSynchronize();
169template <
typename MF,
typename F,
typename FB,
typename P = PostStageNoOp>
170void RK2 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
175 MF dUdt(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
179 fillbndry(1, Uold, time);
180 frhs(1, dUdt, Uold, time,
Real(0.5)*dt);
182 detail::rk_update(Unew, Uold, dUdt, dt);
186 fillbndry(2, Unew, time+dt);
187 frhs(2, dUdt, Unew, time+dt,
Real(0.5)*dt);
190 detail::rk2_update_2(Unew, Uold, dUdt, dt);
207template <
typename MF,
typename F,
typename FB,
typename R,
208 typename P = PostStageNoOp>
209void RK3 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
210 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
215 for (
auto& mf : rkk) {
216 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
221 fillbndry(1, Uold, time);
222 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
224 detail::rk_update(Unew, Uold, rkk[0], dt);
228 fillbndry(2, Unew, time+dt);
229 frhs(2, rkk[1], Unew, time+dt, dt/
Real(6.));
231 detail::rk_update(Unew, Uold, rkk[0], rkk[1],
Real(0.25)*dt);
236 fillbndry(3, Unew, t_half);
237 frhs(3, rkk[2], Unew, t_half, dt*
Real(2./3.));
239 detail::rk3_update_3(Unew, Uold, rkk,
Real(1./6.)*dt);
242 store_crse_data(rkk);
258template <
typename MF,
typename F,
typename FB,
typename R,
259 typename P = PostStageNoOp>
260void RK4 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
261 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
266 for (
auto& mf : rkk) {
267 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
272 fillbndry(1, Uold, time);
273 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
275 detail::rk_update(Unew, Uold, rkk[0],
Real(0.5)*dt);
280 fillbndry(2, Unew, t_half);
281 frhs(2, rkk[1], Unew, t_half, dt/
Real(3.));
283 detail::rk_update(Unew, Uold, rkk[1],
Real(0.5)*dt);
287 fillbndry(3, Unew, t_half);
288 frhs(3, rkk[2], Unew, t_half, dt/
Real(3.));
290 detail::rk_update(Unew, Uold, rkk[2], dt);
294 fillbndry(4, Unew, time+dt);
295 frhs(4, rkk[3], Unew, time+dt, dt/
Real(6.));
297 detail::rk4_update_4(Unew, Uold, rkk,
Real(1./6.)*dt);
300 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:170
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:209
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