Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
AMReX_FEIntegrator.H
Go to the documentation of this file.
1#ifndef AMREX_FE_INTEGRATOR_H
2#define AMREX_FE_INTEGRATOR_H
3#include <AMReX_REAL.H>
4#include <AMReX_Vector.H>
6#include <functional>
7
8namespace amrex {
9
10template<class T>
12{
13private:
15
17
18 // Current (internal) state and time
20 amrex::Real time_current;
21
22 void initialize_stages (const T& S_data, const amrex::Real time)
23 {
24 // Create data for stage RHS
26
27 // Create and initialize data for current state
30
31 // Set the initial time
32 time_current = time;
33 }
34
35public:
37
38 FEIntegrator (const T& S_data, const amrex::Real time = 0.0)
39 {
40 initialize(S_data, time);
41 }
42
43 virtual ~FEIntegrator () {}
44
45 void initialize (const T& S_data, const amrex::Real time = 0.0)
46 {
47 initialize_stages(S_data, time);
48 }
49
50 amrex::Real advance (T& S_old, T& S_new, amrex::Real time, const amrex::Real dt) override
51 {
52 // Assume before step() that S_old is valid data at the current time ("time" argument)
53 // So we initialize S_new by copying the old state.
54 IntegratorOps<T>::Copy(S_new, S_old);
55
56 // F = RHS(S, t)
57 T& F = *F_nodes[0];
58 BaseT::Rhs(F, S_new, time);
59
60 // S_new += timestep * dS/dt
61 IntegratorOps<T>::Saxpy(S_new, dt, F);
62
63 // Call the post step hook
64 BaseT::post_step_action(S_new, time + dt);
65
66 // Return timestep
67 return dt;
68 }
69
70 void evolve (T& S_out, const amrex::Real time_out) override
71 {
72 amrex::Real dt = BaseT::time_step;
73 bool stop = false;
74
75 for (int step_number = 0; step_number < BaseT::max_steps && !stop; ++step_number)
76 {
77 // Adjust step size to reach output time
78 if (time_out - time_current < dt) {
79 dt = time_out - time_current;
80 stop = true;
81 }
82
83 // Call the time integrator step
84 advance(*S_current[0], S_out, time_current, dt);
85
86 // Update current state S_current = S_out
88
89 // Update time
90 time_current += dt;
91
92 if (step_number == BaseT::max_steps - 1) {
93 Error("Did not reach output time in max steps.");
94 }
95 }
96 }
97
98 virtual void time_interpolate (const T& /* S_new */, const T& /* S_old */, amrex::Real /* timestep_fraction */, T& /* data */) override
99 {
100 amrex::Error("Time interpolation not yet supported by the forward euler integrator.");
101 }
102
103 virtual void map_data (std::function<void(T&)> Map) override
104 {
105 for (auto& F : F_nodes) {
106 Map(*F);
107 }
108 }
109
110};
111
112}
113
114#endif
Definition AMReX_FEIntegrator.H:12
void evolve(T &S_out, const amrex::Real time_out) override
Evolve the current (internal) integrator state to time_out.
Definition AMReX_FEIntegrator.H:70
FEIntegrator()
Definition AMReX_FEIntegrator.H:36
amrex::Vector< std::unique_ptr< T > > F_nodes
Definition AMReX_FEIntegrator.H:16
amrex::Vector< std::unique_ptr< T > > S_current
Definition AMReX_FEIntegrator.H:19
amrex::Real advance(T &S_old, T &S_new, amrex::Real time, const amrex::Real dt) override
Take a single time step from (time, S_old) to (time + dt, S_new) with the given step size.
Definition AMReX_FEIntegrator.H:50
FEIntegrator(const T &S_data, const amrex::Real time=0.0)
Definition AMReX_FEIntegrator.H:38
void initialize(const T &S_data, const amrex::Real time=0.0)
Definition AMReX_FEIntegrator.H:45
virtual void time_interpolate(const T &, const T &, amrex::Real, T &) override
Definition AMReX_FEIntegrator.H:98
void initialize_stages(const T &S_data, const amrex::Real time)
Definition AMReX_FEIntegrator.H:22
virtual ~FEIntegrator()
Definition AMReX_FEIntegrator.H:43
virtual void map_data(std::function< void(T &)> Map) override
Definition AMReX_FEIntegrator.H:103
amrex::Real time_current
Definition AMReX_FEIntegrator.H:20
Definition AMReX_IntegratorBase.H:163
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:198
int max_steps
Max number of internal steps before an error is returned (Long)
Definition AMReX_IntegratorBase.H:248
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:168
amrex::Real time_step
Current integrator time step size (Real)
Definition AMReX_IntegratorBase.H:221
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:27
Definition AMReX_Amr.cpp:49
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:224
Definition AMReX_IntegratorBase.H:17