1#ifndef AMREX_INTEGRATOR_BASE_H
2#define AMREX_INTEGRATOR_BASE_H
3#include <AMReX_Config.H>
9#if defined(AMREX_PARTICLES)
20#if defined(AMREX_PARTICLES)
22struct IntegratorOps<T, std::enable_if_t<std::is_base_of_v<amrex::ParticleContainerBase, T> > >
28 V.emplace_back(std::make_unique<T>(Other.Geom(0), Other.ParticleDistributionMap(0), Other.ParticleBoxArray(0)));
29 T& pc = *V[V.size()-1];
36 static void Copy (T& Y,
const T& Other)
39 const bool local =
true;
40 Y.copyParticles(Other, local);
43 static void Saxpy (T& Y,
const amrex::Real a, T& X)
53 auto checkValid = [&]() ->
bool {
54 bool pty_v = pty.isValid();
55 bool ptx_v = ptx.isValid();
57 return pty_v && ptx_v;
60 auto ptIncrement = [&](){ ++pty; ++ptx; };
63#pragma omp parallel if (Gpu::notInLaunchRegion())
65 for (; checkValid(); ptIncrement())
67 const int npy = pty.numParticles();
68 const int npx = ptx.numParticles();
71 ParticleType* psy = &(pty.GetArrayOfStructs()[0]);
72 ParticleType* psx = &(ptx.GetArrayOfStructs()[0]);
74 auto particle_apply_rhs = T::particle_apply_rhs;
77 ParticleType& py = psy[i];
78 const ParticleType& px = psx[i];
79 particle_apply_rhs(py, a, px);
88struct IntegratorOps<T, std::enable_if_t<std::is_same_v<amrex::Vector<amrex::MultiFab>, T> > >
94 V.emplace_back(std::make_unique<T>());
95 for (
auto const& other_mf : Other) {
97 V.back()->push_back(
amrex::MultiFab(other_mf.boxArray(), other_mf.DistributionMap(), other_mf.nComp(), nGrow));
104 const int size = Y.size();
105 bool specify_components = !scomp.empty() && ncomp.
size() == scomp.size();
106 for (
int i = 0; i < size; ++i) {
108 const int iscomp = specify_components ? scomp[i] : 0;
109 const int incomp = specify_components ? ncomp[i] : Other[i].nComp();
119 const int size = Y.size();
120 bool specify_components = !scomp.empty() && ncomp.
size() == scomp.size();
121 for (
int i = 0; i < size; ++i) {
123 const int iscomp = specify_components ? scomp[i] : 0;
124 const int incomp = specify_components ? ncomp[i] : X[i].nComp();
134struct IntegratorOps<T, std::enable_if_t<std::is_same_v<amrex::MultiFab, T> > >
141 V.emplace_back(std::make_unique<T>(Other.boxArray(), Other.DistributionMap(), Other.nComp(), nGrow));
144 static void Copy (T& Y,
const T& Other,
const int scomp=0,
const int ncomp=-1,
bool Grow =
true)
148 const int mf_ncomp = ncomp > 0 ? ncomp : Other.nComp();
152 static void Saxpy (T& Y,
const amrex::Real a,
const T& X,
const int scomp=0,
const int ncomp=-1,
bool Grow =
false)
156 const int mf_ncomp = ncomp > 0 ? ncomp : X.nComp();
182 std::function<void(T& rhs, T& state,
const amrex::Real time)>
Rhs;
188 std::function<void(T& rhs, T& state,
const amrex::Real time)>
RhsIm;
194 std::function<void(T& rhs, T& state,
const amrex::Real time)>
RhsEx;
200 std::function<void(T& rhs, T& state,
const amrex::Real time)>
RhsFast;
298 void set_rhs (std::function<
void(T&, T&,
const amrex::Real)>
F)
304 std::function<
void(T&, T&,
const amrex::Real)> Fe)
383 virtual amrex::Real
advance (T& S_old, T& S_new, amrex::Real time,
389 virtual void evolve (T& S_out,
const amrex::Real time_out) = 0;
392 amrex::Real timestep_fraction, T& data) = 0;
394 virtual void map_data (std::function<
void(T&)> Map) = 0;
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
Definition AMReX_IntegratorBase.H:164
void set_post_step_action(std::function< void(T &, amrex::Real)> A)
Definition AMReX_IntegratorBase.H:320
void set_rhs(std::function< void(T &, T &, const amrex::Real)> F)
Definition AMReX_IntegratorBase.H:298
void set_adaptive_step()
Definition AMReX_IntegratorBase.H:346
bool use_adaptive_fast_time_step
Flag to enable/disable adaptive time stepping at the fast time scale in multirate methods (bool)
Definition AMReX_IntegratorBase.H:246
amrex::Real fast_rel_tol
Relative tolerance for adaptive time stepping at the fast time scale (Real)
Definition AMReX_IntegratorBase.H:278
amrex::Real rel_tol
Relative tolerance for adaptive time stepping (Real)
Definition AMReX_IntegratorBase.H:267
void initialize_base_parameters()
Definition AMReX_IntegratorBase.H:166
void set_fast_tolerances(amrex::Real rtol, amrex::Real atol)
Definition AMReX_IntegratorBase.H:373
void set_post_fast_stage_action(std::function< void(T &, amrex::Real)> A)
Definition AMReX_IntegratorBase.H:325
std::function< void(T &, amrex::Real)> post_fast_stage_action
The post_stage_action function is called by the integrator on the computed stage just after it is com...
Definition AMReX_IntegratorBase.H:218
virtual void time_interpolate(const T &S_new, const T &S_old, amrex::Real timestep_fraction, T &data)=0
void set_adaptive_fast_step()
Definition AMReX_IntegratorBase.H:357
virtual amrex::Real advance(T &S_old, T &S_new, amrex::Real time, amrex::Real dt)=0
Take a single time step from (time, S_old) to (time + dt, S_new) with the given step size.
amrex::Real get_time_step()
Definition AMReX_IntegratorBase.H:335
amrex::Real fast_abs_tol
Absolute tolerance for adaptive time stepping at the fast time scale (Real)
Definition AMReX_IntegratorBase.H:284
amrex::Real fast_time_step
Current integrator fast time scale time step size with multirate methods (Real)
Definition AMReX_IntegratorBase.H:252
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsEx
RhsEx is the explicit right-hand-side function an ImEx integrator will use.
Definition AMReX_IntegratorBase.H:194
std::function< void(T &, amrex::Real)> post_step_action
The post_step_action function is called by the integrator on the computed state just after it is comp...
Definition AMReX_IntegratorBase.H:212
void set_post_stage_action(std::function< void(T &, amrex::Real)> A)
Definition AMReX_IntegratorBase.H:315
IntegratorBase()
Definition AMReX_IntegratorBase.H:288
std::function< void(T &, amrex::Real)> post_fast_step_action
The post_step_action function is called by the integrator on the computed state just after it is comp...
Definition AMReX_IntegratorBase.H:224
int max_steps
Max number of internal steps before an error is returned (Long)
Definition AMReX_IntegratorBase.H:262
virtual void evolve(T &S_out, const amrex::Real time_out)=0
Evolve the current (internal) integrator state to time_out.
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsIm
RhsIm is the implicit right-hand-side function an ImEx integrator will use.
Definition AMReX_IntegratorBase.H:188
virtual void map_data(std::function< void(T &)> Map)=0
void set_post_fast_step_action(std::function< void(T &, amrex::Real)> A)
Definition AMReX_IntegratorBase.H:330
void set_imex_rhs(std::function< void(T &, T &, const amrex::Real)> Fi, std::function< void(T &, T &, const amrex::Real)> Fe)
Definition AMReX_IntegratorBase.H:303
void set_max_steps(int steps)
Definition AMReX_IntegratorBase.H:362
void set_time_step(amrex::Real dt)
Definition AMReX_IntegratorBase.H:340
void set_tolerances(amrex::Real rtol, amrex::Real atol)
Definition AMReX_IntegratorBase.H:367
std::function< void(T &rhs, T &state, const amrex::Real time)> Rhs
Rhs is the right-hand-side function the integrator will use.
Definition AMReX_IntegratorBase.H:182
virtual ~IntegratorBase()=default
amrex::Long num_steps
Number of integrator time steps (Long)
Definition AMReX_IntegratorBase.H:257
bool use_adaptive_time_step
Flag to enable/disable adaptive time stepping in single rate methods or at the slow time scale in mul...
Definition AMReX_IntegratorBase.H:230
void set_fast_rhs(std::function< void(T &, T &, const amrex::Real)> F)
Definition AMReX_IntegratorBase.H:310
std::function< void(T &, amrex::Real)> post_stage_action
The post_stage_action function is called by the integrator on the computed stage just after it is com...
Definition AMReX_IntegratorBase.H:206
void set_fast_time_step(amrex::Real dt)
Definition AMReX_IntegratorBase.H:351
IntegratorBase(const T &)
Definition AMReX_IntegratorBase.H:292
amrex::Real time_step
Current integrator time step size (Real)
Definition AMReX_IntegratorBase.H:235
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsFast
RhsFast is the fast timescale right-hand-side function a multirate integrator will use.
Definition AMReX_IntegratorBase.H:200
amrex::Real previous_time_step
Step size of the last completed step (Real)
Definition AMReX_IntegratorBase.H:240
amrex::Real abs_tol
Absolute tolerance for adaptive time stepping (Real)
Definition AMReX_IntegratorBase.H:272
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
static void Saxpy(MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
dst += a*src
Definition AMReX_MultiFab.cpp:292
static void Copy(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Copy from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying Bo...
Definition AMReX_MultiFab.cpp:193
Definition AMReX_ParIter.H:115
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:343
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1393
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
Definition AMReX_Amr.cpp:49
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:191
void Copy(FabArray< DFAB > &dst, FabArray< SFAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:180
IntVectND< 3 > IntVect
Definition AMReX_BaseFwd.H:30
static void CreateLike(amrex::Vector< std::unique_ptr< T > > &V, const T &Other)
Definition AMReX_IntegratorBase.H:25
static void Saxpy(T &Y, const amrex::Real a, T &X)
Definition AMReX_IntegratorBase.H:43
static void Copy(T &Y, const T &Other)
Definition AMReX_IntegratorBase.H:36
static void CreateLike(amrex::Vector< std::unique_ptr< T > > &V, const T &Other, bool Grow=false)
Definition AMReX_IntegratorBase.H:137
static void Saxpy(T &Y, const amrex::Real a, const T &X, const int scomp=0, const int ncomp=-1, bool Grow=false)
Definition AMReX_IntegratorBase.H:152
static void Copy(T &Y, const T &Other, const int scomp=0, const int ncomp=-1, bool Grow=true)
Definition AMReX_IntegratorBase.H:144
static void Saxpy(T &Y, const amrex::Real a, const T &X, const Vector< int > &scomp={}, const Vector< int > &ncomp={}, bool Grow=false)
Definition AMReX_IntegratorBase.H:116
static void CreateLike(amrex::Vector< std::unique_ptr< T > > &V, const T &Other, bool Grow=false)
Definition AMReX_IntegratorBase.H:91
static void Copy(T &Y, const T &Other, const Vector< int > &scomp={}, const Vector< int > &ncomp={}, bool Grow=true)
Definition AMReX_IntegratorBase.H:101
Definition AMReX_IntegratorBase.H:18
The struct used to store particles.
Definition AMReX_Particle.H:402