Block-Structured AMR Software Framework
 
Loading...
Searching...
No Matches
AMReX_RungeKutta.H
Go to the documentation of this file.
1#ifndef AMREX_RUNGE_KUTTA_H_
2#define AMREX_RUNGE_KUTTA_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_FabArray.H>
6
50
52 template <typename MF>
53 std::enable_if_t<IsFabArray<MF>::value> operator() (int, MF&) const {}
54};
55
57namespace detail {
59template <typename MF>
60void rk_update (MF& Unew, MF const& Uold, MF const& dUdt, Real dt)
61{
62 auto const& snew = Unew.arrays();
63 auto const& sold = Uold.const_arrays();
64 auto const& sdot = dUdt.const_arrays();
65 amrex::ParallelFor(Unew, IntVect(0), Unew.nComp(), [=] AMREX_GPU_DEVICE
66 (int bi, int i, int j, int k, int n) noexcept
67 {
68 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*sdot[bi](i,j,k,n);
69 });
70 Gpu::streamSynchronize();
71}
72
74template <typename MF>
75void rk_update (MF& Unew, MF const& Uold, MF const& dUdt1, MF const& dUdt2, Real dt)
76{
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();
81 amrex::ParallelFor(Unew, IntVect(0), Unew.nComp(), [=] AMREX_GPU_DEVICE
82 (int bi, int i, int j, int k, int n) noexcept
83 {
84 snew[bi](i,j,k,n) = sold[bi](i,j,k,n) + dt*(sdot1[bi](i,j,k,n) +
85 sdot2[bi](i,j,k,n));
86 });
87 Gpu::streamSynchronize();
88}
89
91template <typename MF>
92void rk2_update_2 (MF& Unew, MF const& Uold, MF const& dUdt, Real dt)
93{
94 auto const& snew = Unew.arrays();
95 auto const& sold = Uold.const_arrays();
96 auto const& sdot = dUdt.const_arrays();
97 amrex::ParallelFor(Unew, IntVect(0), Unew.nComp(), [=] AMREX_GPU_DEVICE
98 (int bi, int i, int j, int k, int n) noexcept
99 {
100 snew[bi](i,j,k,n) = Real(0.5)*(snew[bi](i,j,k,n) +
101 sold[bi](i,j,k,n) +
102 sdot[bi](i,j,k,n) * dt);
103 });
104 Gpu::streamSynchronize();
105}
106
108template <typename MF>
109void rk3_update_3 (MF& Unew, MF const& Uold, Array<MF,3> const& rkk, Real dt6)
110{
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();
116 amrex::ParallelFor(Unew, IntVect(0), Unew.nComp(), [=] AMREX_GPU_DEVICE
117 (int bi, int i, int j, int k, int n) noexcept
118 {
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));
122 });
123 Gpu::streamSynchronize();
124}
125
127template <typename MF>
128void rk4_update_4 (MF& Unew, MF const& Uold, Array<MF,4> const& rkk, Real dt6)
129{
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();
136 amrex::ParallelFor(Unew, IntVect(0), Unew.nComp(), [=] AMREX_GPU_DEVICE
137 (int bi, int i, int j, int k, int n) noexcept
138 {
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)));
142 });
143 Gpu::streamSynchronize();
144}
145}
147
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,
161 P const& post_stage = PostStageNoOp())
162{
163 BL_PROFILE("RungeKutta2");
164
165 MF dUdt(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
166 MFInfo(), Unew.Factory());
167
168 // RK2 stage 1
169 fillbndry(1, Uold, time);
170 frhs(1, dUdt, Uold, time, Real(0.5)*dt);
171 // Unew = Uold + dt * dUdt
172 detail::rk_update(Unew, Uold, dUdt, dt);
173 post_stage(1, Unew);
174
175 // RK2 stage 2
176 fillbndry(2, Unew, time+dt);
177 frhs(2, dUdt, Unew, time, Real(0.5)*dt);
178 // Unew = (Uold+Unew)/2 + dUdt_2 * dt/2,
179 // which is Unew = Uold + dt/2 * (dUdt_1 + dUdt_2)
180 detail::rk2_update_2(Unew, Uold, dUdt, dt);
181 post_stage(2, Unew);
182}
183
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())
200{
201 BL_PROFILE("RungeKutta3");
202
203 Array<MF,3> rkk;
204 for (auto& mf : rkk) {
205 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
206 MFInfo(), Unew.Factory());
207 }
208
209 // RK3 stage 1
210 fillbndry(1, Uold, time);
211 frhs(1, rkk[0], Uold, time, dt/Real(6.));
212 // Unew = Uold + k1 * dt
213 detail::rk_update(Unew, Uold, rkk[0], dt);
214 post_stage(1, Unew);
215
216 // RK3 stage 2
217 fillbndry(2, Unew, time+dt);
218 frhs(2, rkk[1], Unew, time+dt, dt/Real(6.));
219 // Unew = Uold + (k1+k2) * dt/4
220 detail::rk_update(Unew, Uold, rkk[0], rkk[1], Real(0.25)*dt);
221 post_stage(2, Unew);
222
223 // RK3 stage 3
224 Real t_half = time + Real(0.5)*dt;
225 fillbndry(3, Unew, t_half);
226 frhs(3, rkk[2], Unew, t_half, dt*Real(2./3.));
227 // Unew = Uold + (k1/6 + k2/6 + k3*(2/3)) * dt
228 detail::rk3_update_3(Unew, Uold, rkk, Real(1./6.)*dt);
229 post_stage(3, Unew);
230
231 store_crse_data(rkk);
232}
233
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())
250{
251 BL_PROFILE("RungeKutta4");
252
253 Array<MF,4> rkk;
254 for (auto& mf : rkk) {
255 mf.define(Unew.boxArray(), Unew.DistributionMap(), Unew.nComp(), 0,
256 MFInfo(), Unew.Factory());
257 }
258
259 // RK4 stage 1
260 fillbndry(1, Uold, time);
261 frhs(1, rkk[0], Uold, time, dt/Real(6.));
262 // Unew = Uold + k1 * dt/2
263 detail::rk_update(Unew, Uold, rkk[0], Real(0.5)*dt);
264 post_stage(1, Unew);
265
266 // RK4 stage 2
267 Real t_half = time + Real(0.5)*dt;
268 fillbndry(2, Unew, t_half);
269 frhs(2, rkk[1], Unew, t_half, dt/Real(3.));
270 // Unew = Uold + k2 * dt/2
271 detail::rk_update(Unew, Uold, rkk[1], Real(0.5)*dt);
272 post_stage(2, Unew);
273
274 // RK4 stage 3
275 fillbndry(3, Unew, t_half);
276 frhs(3, rkk[2], Unew, t_half, dt/Real(3.));
277 // Unew = Uold + k3 * dt;
278 detail::rk_update(Unew, Uold, rkk[2], dt);
279 post_stage(3, Unew);
280
281 // RK4 stage 4
282 fillbndry(4, Unew, time+dt);
283 frhs(4, rkk[3], Unew, time+dt, dt/Real(6.));
284 // Unew = Uold + (k1/6 + k2/3 + k3/3 + k4/6) * dt
285 detail::rk4_update_4(Unew, Uold, rkk, Real(1./6.)*dt);
286 post_stage(4, Unew);
287
288 store_crse_data(rkk);
289}
290
291}
292
293#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:25
Functions for Runge-Kutta methods.
Definition AMReX_RungeKutta.H:49
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
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_RungeKutta.H:51
std::enable_if_t< IsFabArray< MF >::value > operator()(int, MF &) const
Definition AMReX_RungeKutta.H:53