PFASST++
pfasst::examples::scalar::ScalarSweeper< time > Class Template Reference

Sweeper for scalar test equation. More...

#include <scalar_sweeper.hpp>

+ Inheritance diagram for pfasst::examples::scalar::ScalarSweeper< time >:
+ Collaboration diagram for pfasst::examples::scalar::ScalarSweeper< time >:

Public Member Functions

 ScalarSweeper (const complex< double > &lambda, const complex< double > &u0)
 generic constructor; initialize all function call counters with zero. More...
 
virtual ~ScalarSweeper ()
 upon destruction, report final error and number of function calls More...
 
void echo_error (time t)
 compute error between last state and exact solution at time tand print it to cout More...
 
double get_errors ()
 returns error, but does not update it! More...
 
void post_predict () override
 post prediction step, update of error. More...
 
void post_sweep () override
 post sweep, update error. More...
 
void exact (complex_vector_type &q, time t)
 computes the exact solution \( u_0 \exp \left( \lambda*t \right) \) at a given time t. More...
 
void exact (shared_ptr< encap_type > q_encap, time t)
 
void f_expl_eval (shared_ptr< encap_type > f_encap, shared_ptr< encap_type > q_encap, time t) override
 evaluate the explicit part of the right hand side: Multiply with \( \text{imag}(\lambda) \) More...
 
void f_impl_eval (shared_ptr< encap_type > f_encap, shared_ptr< encap_type > q_encap, time t) override
 evaluate the implicit part of the right hand side: Multiply with \( \text{real}(\lambda) \) More...
 
void impl_solve (shared_ptr< encap_type > f_encap, shared_ptr< encap_type > q_encap, time t, time dt, shared_ptr< encap_type > rhs_encap) override
 for given \( b \), solve \( \left( \mathbb{I}_d - \Delta t \text{real}(\lambda) \right) u = b \) for \( u \) and set f_encap to \( \text{real}(\lambda) u \) More...
 
- Public Member Functions inherited from pfasst::encap::IMEXSweeper< time >
 IMEXSweeper ()=default
 
virtual ~IMEXSweeper ()=default
 
virtual void setup (bool coarse) override
 Setup (allocate etc) the sweeper. More...
 
virtual void predict (bool initial) override
 Compute low-order provisional solution. More...
 
virtual void sweep () override
 Perform one SDC sweep/iteration. More...
 
virtual void advance () override
 Advance the end solution to start solution. More...
 
virtual void reevaluate (bool initial_only) override
 Re-evaluate function values. More...
 
virtual void integrate (time dt, vector< shared_ptr< Encapsulation< time >>> dst) const override
 Integrates values of right hand side at all time nodes \( t \in [0,M-1] \) simultaneously. More...
 
virtual void residual (time dt, vector< shared_ptr< Encapsulation< time >>> dst) const override
 Compute residual at each SDC node (including FAS corrections). More...
 
virtual void f_expl_eval (shared_ptr< Encapsulation< time >> f_expl_encap, shared_ptr< Encapsulation< time >> u_encap, time t)
 Evaluate the explicit part of the ODE. More...
 
virtual void f_impl_eval (shared_ptr< Encapsulation< time >> f_impl_encap, shared_ptr< Encapsulation< time >> u_encap, time t)
 Evaluate the implicit part of the ODE. More...
 
virtual void impl_solve (shared_ptr< Encapsulation< time >> f_encap, shared_ptr< Encapsulation< time >> u_encap, time t, time dt, shared_ptr< Encapsulation< time >> rhs_encap)
 Solve \( U - \Delta t F_{\rm impl}(U) = RHS \) for \( U \). More...
 
- Public Member Functions inherited from pfasst::encap::EncapSweeper< time >
 EncapSweeper ()
 
virtual shared_ptr< Encapsulation< time > > get_state (size_t m) const
 Retrieve solution values of current iteration at time node index m. More...
 
virtual shared_ptr< Encapsulation< time > > get_tau (size_t m) const
 Retrieve FAS correction of current iteration at time node index m. More...
 
virtual shared_ptr< Encapsulation< time > > get_saved_state (size_t m) const
 Retrieve solution values of previous iteration at time node index m. More...
 
virtual void set_options () override
 Set options from command line etc. More...
 
virtual void spread () override
 Initialize solution values at all time nodes with meaningful values. More...
 
virtual void save (bool initial_only) override
 Save states (and/or function values) at all nodes. More...
 
virtual void set_quadrature (shared_ptr< IQuadrature< time >> quadrature)
 
virtual shared_ptr< const IQuadrature< time > > get_quadrature () const
 
virtual const vector< time > get_nodes () const
 
virtual void set_factory (shared_ptr< EncapFactory< time >> factory)
 
virtual shared_ptr< EncapFactory< time > > get_factory () const
 
virtual shared_ptr< Encapsulation< time > > get_start_state () const
 
virtual shared_ptr< Encapsulation< time > > get_end_state ()
 
void set_residual_tolerances (time abs_residual_tol, time rel_residual_tol, int order=0)
 Set residual tolerances for convergence checking. More...
 
virtual bool converged () override
 Return convergence status. More...
 
virtual void post (ICommunicator *comm, int tag) override
 
virtual void send (ICommunicator *comm, int tag, bool blocking) override
 
virtual void recv (ICommunicator *comm, int tag, bool blocking) override
 
virtual void broadcast (ICommunicator *comm) override
 
- Public Member Functions inherited from pfasst::ISweeper< time >
 ISweeper ()
 
virtual ~ISweeper ()
 
virtual void set_controller (Controller< time > *ctrl)
 Set the sweepers controller. More...
 
virtual Controller< time > * get_controller ()
 Accessor to the controller managing this sweeper. More...
 
virtual void post_step ()
 Hook automatically run after each completed time step. More...
 

Private Types

typedef encap::Encapsulation< time > encap_type
 
typedef encap::VectorEncapsulation< complex< double > > complex_vector_type
 define a type for a complex PFASST vector encapsulation. More...
 

Private Attributes

complex< double > lambda
 parameter lambda and initial value \( u_0 \) More...
 
complex< double > u0
 
const complex< double > i_complex = complex<double>(0, 1)
 the complex unit \( i = \sqrt{-1} \) More...
 
double error
 error at the final time. For the scalar example, an analytical solution is known. More...
 
size_t n_f_expl_eval
 counters for how often f_expl_eval, f_impl_eval and impl_solve are called. More...
 
size_t n_f_impl_eval
 
size_t n_impl_solve
 

Additional Inherited Members

- Protected Member Functions inherited from pfasst::encap::IMEXSweeper< time >
virtual void integrate_end_state (time dt)
 Set end state to \( U_0 + \int F_{expl} + F_{expl} \). More...
 
- Protected Attributes inherited from pfasst::encap::IMEXSweeper< time >
vector< shared_ptr< Encapsulation< time > > > s_integrals
 Node-to-node integrals of \( F(t,u) \) at all time nodes of the current iteration. More...
 
vector< shared_ptr< Encapsulation< time > > > fs_expl
 Values of the explicit part of the right hand side \( F_{expl}(t,u) \) at all time nodes of the current iteration. More...
 
shared_ptr< Encapsulation< time > > fs_expl_start
 
vector< shared_ptr< Encapsulation< time > > > fs_impl
 Values of the implicit part of the right hand side \( F_{impl}(t,u) \) at all time nodes of the current iteration. More...
 
- Protected Attributes inherited from pfasst::encap::EncapSweeper< time >
string FORMAT_STR
 
shared_ptr< IQuadrature< time > > quadrature
 
shared_ptr< EncapFactory< time > > factory
 Encapsulation data structure factory. More...
 
shared_ptr< Encapsulation< time > > start_state
 Separate start state, i.e. initial condition for the sweeper's current time step. More...
 
shared_ptr< Encapsulation< time > > end_state
 Current solution at \( T_{end} \). More...
 
vector< shared_ptr< Encapsulation< time > > > residuals
 Place for the residuals at the different time nodes. More...
 
vector< shared_ptr< Encapsulation< time > > > state
 Solution values \( U \) at all time nodes of the current iteration. More...
 
vector< shared_ptr< Encapsulation< time > > > saved_state
 Solution values \( U \) at all time nodes of the previous iteration. More...
 
vector< shared_ptr< Encapsulation< time > > > fas_corrections
 FAS corrections \( \tau \) at all time nodes of the current iteration. More...
 
int residual_norm_order
 
time abs_residual_tol
 Tolerance for absolute residual. More...
 
time rel_residual_tol
 Tolerance for relative residual. More...
 
- Protected Attributes inherited from pfasst::ISweeper< time >
Controller< time > * controller
 Backreference to the controller managing the sweeper instance. More...
 

Detailed Description

template<typename time = pfasst::time_precision>
class pfasst::examples::scalar::ScalarSweeper< time >

Sweeper for scalar test equation.

\[ u' = \lambda*u \quad\text{ , } u(0) = u_0 \] with complex lambda using an IMEX scheme.

Definition at line 51 of file scalar_sweeper.hpp.

Member Typedef Documentation

template<typename time = pfasst::time_precision>
typedef encap::VectorEncapsulation<complex<double> > pfasst::examples::scalar::ScalarSweeper< time >::complex_vector_type
private

define a type for a complex PFASST vector encapsulation.

Definition at line 58 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
typedef encap::Encapsulation<time> pfasst::examples::scalar::ScalarSweeper< time >::encap_type
private

Definition at line 55 of file scalar_sweeper.hpp.

Constructor & Destructor Documentation

template<typename time = pfasst::time_precision>
pfasst::examples::scalar::ScalarSweeper< time >::ScalarSweeper ( const complex< double > &  lambda,
const complex< double > &  u0 
)
inline

generic constructor; initialize all function call counters with zero.

Parameters
[in]lambdacoefficient in test equation
[in]u0initial value at \( t=0 \)

Definition at line 79 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
virtual pfasst::examples::scalar::ScalarSweeper< time >::~ScalarSweeper ( )
inlinevirtual

upon destruction, report final error and number of function calls

Definition at line 91 of file scalar_sweeper.hpp.

References ML_LOG.

Member Function Documentation

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::echo_error ( time  t)
inline

compute error between last state and exact solution at time tand print it to cout

Parameters
[in]tTime

Definition at line 104 of file scalar_sweeper.hpp.

References ML_LOG.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::exact ( complex_vector_type q,
time  t 
)
inline

computes the exact solution \( u_0 \exp \left( \lambda*t \right) \) at a given time t.

Parameters
[in]tTime

Definition at line 149 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::exact ( shared_ptr< encap_type q_encap,
time  t 
)
inline

Definition at line 154 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::f_expl_eval ( shared_ptr< encap_type f_encap,
shared_ptr< encap_type q_encap,
time  t 
)
inlineoverride

evaluate the explicit part of the right hand side: Multiply with \( \text{imag}(\lambda) \)

Definition at line 163 of file scalar_sweeper.hpp.

References UNUSED.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::f_impl_eval ( shared_ptr< encap_type f_encap,
shared_ptr< encap_type q_encap,
time  t 
)
inlineoverride

evaluate the implicit part of the right hand side: Multiply with \( \text{real}(\lambda) \)

Definition at line 179 of file scalar_sweeper.hpp.

References UNUSED.

template<typename time = pfasst::time_precision>
double pfasst::examples::scalar::ScalarSweeper< time >::get_errors ( )
inline

returns error, but does not update it!

Definition at line 119 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::impl_solve ( shared_ptr< encap_type f_encap,
shared_ptr< encap_type q_encap,
time  t,
time  dt,
shared_ptr< encap_type rhs_encap 
)
inlineoverride

for given \( b \), solve \( \left( \mathbb{I}_d - \Delta t \text{real}(\lambda) \right) u = b \) for \( u \) and set f_encap to \( \text{real}(\lambda) u \)

Definition at line 197 of file scalar_sweeper.hpp.

References UNUSED.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::post_predict ( )
inlineoverridevirtual

post prediction step, update of error.

Reimplemented from pfasst::ISweeper< time >.

Definition at line 127 of file scalar_sweeper.hpp.

References plot::dt, and plot::t.

template<typename time = pfasst::time_precision>
void pfasst::examples::scalar::ScalarSweeper< time >::post_sweep ( )
inlineoverridevirtual

post sweep, update error.

Reimplemented from pfasst::ISweeper< time >.

Definition at line 137 of file scalar_sweeper.hpp.

References plot::dt, and plot::t.

Member Data Documentation

template<typename time = pfasst::time_precision>
double pfasst::examples::scalar::ScalarSweeper< time >::error
private

error at the final time. For the scalar example, an analytical solution is known.

Definition at line 67 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
const complex<double> pfasst::examples::scalar::ScalarSweeper< time >::i_complex = complex<double>(0, 1)
private

the complex unit \( i = \sqrt{-1} \)

Definition at line 64 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
complex<double> pfasst::examples::scalar::ScalarSweeper< time >::lambda
private

parameter lambda and initial value \( u_0 \)

Definition at line 61 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::scalar::ScalarSweeper< time >::n_f_expl_eval
private

counters for how often f_expl_eval, f_impl_eval and impl_solve are called.

Definition at line 70 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::scalar::ScalarSweeper< time >::n_f_impl_eval
private

Definition at line 70 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::scalar::ScalarSweeper< time >::n_impl_solve
private

Definition at line 70 of file scalar_sweeper.hpp.

template<typename time = pfasst::time_precision>
complex<double> pfasst::examples::scalar::ScalarSweeper< time >::u0
private

Definition at line 61 of file scalar_sweeper.hpp.


The documentation for this class was generated from the following file: