53 template <FabArrayType MF>
61void rk_update (MF& Unew, MF
const& Uold, MF
const& dUdt,
Real dt)
63 auto const& snew = Unew.arrays();
64 auto const& sold = Uold.const_arrays();
65 auto const& sdot = dUdt.const_arrays();
67 (
int bi,
int i,
int j,
int k,
int n)
noexcept
69 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*sdot[bi](i,j,k,n);
71 if (!Gpu::inNoSyncRegion()) {
72 Gpu::streamSynchronize();
78void rk_update (MF& Unew, MF
const& Uold, MF
const& dUdt1, MF
const& dUdt2, Real dt)
80 auto const& snew = Unew.arrays();
81 auto const& sold = Uold.const_arrays();
82 auto const& sdot1 = dUdt1.const_arrays();
83 auto const& sdot2 = dUdt2.const_arrays();
85 (
int bi,
int i,
int j,
int k,
int n)
noexcept
87 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*(sdot1[bi](i,j,k,n) +
90 if (!Gpu::inNoSyncRegion()) {
91 Gpu::streamSynchronize();
97void rk2_update_2 (MF& Unew, MF
const& Uold, MF
const& dUdt, Real dt)
99 auto const& snew = Unew.arrays();
100 auto const& sold = Uold.const_arrays();
101 auto const& sdot = dUdt.const_arrays();
103 (
int bi,
int i,
int j,
int k,
int n)
noexcept
105 snew[bi](i,j,k,n) =
Real(0.5)*(snew[bi](i,j,k,n) +
107 sdot[bi](i,j,k,n) * dt);
109 if (!Gpu::inNoSyncRegion()) {
110 Gpu::streamSynchronize();
115template <
typename MF>
116void rk3_update_3 (MF& Unew, MF
const& Uold, Array<MF,3>
const& rkk, Real dt6)
118 auto const& snew = Unew.arrays();
119 auto const& sold = Uold.const_arrays();
120 auto const& k1 = rkk[0].const_arrays();
121 auto const& k2 = rkk[1].const_arrays();
122 auto const& k3 = rkk[2].const_arrays();
124 (
int bi,
int i,
int j,
int k,
int n)
noexcept
126 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
127 + dt6 * (k1[bi](i,j,k,n) + k2[bi](i,j,k,n)
128 +
Real(4.) * k3[bi](i,j,k,n));
130 if (!Gpu::inNoSyncRegion()) {
131 Gpu::streamSynchronize();
136template <
typename MF>
137void rk4_update_4 (MF& Unew, MF
const& Uold, Array<MF,4>
const& rkk, Real dt6)
139 auto const& snew = Unew.arrays();
140 auto const& sold = Uold.const_arrays();
141 auto const& k1 = rkk[0].const_arrays();
142 auto const& k2 = rkk[1].const_arrays();
143 auto const& k3 = rkk[2].const_arrays();
144 auto const& k4 = rkk[3].const_arrays();
146 (
int bi,
int i,
int j,
int k,
int n)
noexcept
148 snew[bi](i,j,k,n) = sold[bi](i,j,k,n)
149 + dt6 * ( k1[bi](i,j,k,n) + k4[bi](i,j,k,n)
150 +
Real(2.)*(k2[bi](i,j,k,n) + k3[bi](i,j,k,n)));
152 if (!Gpu::inNoSyncRegion()) {
153 Gpu::streamSynchronize();
170template <
typename MF,
typename F,
typename FB,
typename P = PostStageNoOp>
171void RK2 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
176 MF dUdt(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
180 fillbndry(1, Uold, time);
181 frhs(1, dUdt, Uold, time,
Real(0.5)*dt);
183 detail::rk_update(Unew, Uold, dUdt, dt);
187 fillbndry(2, Unew, time+dt);
188 frhs(2, dUdt, Unew, time+dt,
Real(0.5)*dt);
191 detail::rk2_update_2(Unew, Uold, dUdt, dt);
208template <
typename MF,
typename F,
typename FB,
typename R,
209 typename P = PostStageNoOp>
210void RK3 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
211 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
216 for (
auto& mf : rkk) {
217 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
222 fillbndry(1, Uold, time);
223 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
225 detail::rk_update(Unew, Uold, rkk[0], dt);
229 fillbndry(2, Unew, time+dt);
230 frhs(2, rkk[1], Unew, time+dt, dt/
Real(6.));
232 detail::rk_update(Unew, Uold, rkk[0], rkk[1],
Real(0.25)*dt);
237 fillbndry(3, Unew, t_half);
238 frhs(3, rkk[2], Unew, t_half, dt*
Real(2./3.));
240 detail::rk3_update_3(Unew, Uold, rkk,
Real(1./6.)*dt);
243 store_crse_data(rkk);
259template <
typename MF,
typename F,
typename FB,
typename R,
260 typename P = PostStageNoOp>
261void RK4 (MF& Uold, MF& Unew,
Real time,
Real dt,
F const& frhs, FB
const& fillbndry,
262 R
const& store_crse_data, P
const& post_stage =
PostStageNoOp())
267 for (
auto& mf : rkk) {
268 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
273 fillbndry(1, Uold, time);
274 frhs(1, rkk[0], Uold, time, dt/
Real(6.));
276 detail::rk_update(Unew, Uold, rkk[0],
Real(0.5)*dt);
281 fillbndry(2, Unew, t_half);
282 frhs(2, rkk[1], Unew, t_half, dt/
Real(3.));
284 detail::rk_update(Unew, Uold, rkk[1],
Real(0.5)*dt);
288 fillbndry(3, Unew, t_half);
289 frhs(3, rkk[2], Unew, t_half, dt/
Real(3.));
291 detail::rk_update(Unew, Uold, rkk[2], dt);
295 fillbndry(4, Unew, time+dt);
296 frhs(4, rkk[3], Unew, time+dt, dt/
Real(6.));
298 detail::rk4_update_4(Unew, Uold, rkk,
Real(1./6.)*dt);
301 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:171
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:210
Definition AMReX_RungeKutta.H:52
void operator()(int, MF &) const
Definition AMReX_RungeKutta.H:54