PFASST++
pfasst::examples::vdp::VdpSweeper< time > Class Template Reference

Sweeper for the van der Pol oscillator. More...

#include <vdp_sweeper.hpp>

+ Inheritance diagram for pfasst::examples::vdp::VdpSweeper< time >:
+ Collaboration diagram for pfasst::examples::vdp::VdpSweeper< time >:

Public Member Functions

 VdpSweeper (double nu, double x0, double y0)
 generic constructor; initialize all function call counters with zero. More...
 
virtual ~VdpSweeper ()
 upon destruction, report final error and number of function calls and close output file. More...
 
void echo_error (time t)
 computes and prints out error, which is meaningful only for \( \nu=0 \). More...
 
double get_errors ()
 returns error, but does not update it! More...
 
void post_predict () override
 post prediction step, update error. More...
 
void post_sweep () override
 post sweep, update error. More...
 
void exact (real_vector_type &q, time t)
 if \( \nu=0 \), this computes the exact solution at time t. More...
 
void exact (shared_ptr< encap_type > q_encap, time t)
 
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. 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( u - \Delta t f(u) \right) u = b \) for \( u \) and set f_encap to \( f(u) \). More...
 
- Public Member Functions inherited from pfasst::encap::ImplicitSweeper< time >
 ImplicitSweeper ()=default
 
virtual ~ImplicitSweeper ()=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
 
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 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 void residual (time dt, vector< shared_ptr< Encapsulation< time >>> dst) const
 Compute residual at each SDC node (including FAS corrections). 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< double > real_vector_type
 Define a type for a complex PFASST vector encapsulation. More...
 

Private Attributes

double nu
 Parameter \( \nu \) in the van-der-Pol oscillator. More...
 
double x0
 Starting values. More...
 
double y0
 
size_t newton_maxit
 Maximum number of iterations for the nonlinear Newton solver. More...
 
double newton_tol
 Tolerance for the nonlinear Newton iteration. More...
 
size_t n_f_impl_eval
 Counters for how often f_impl_eval and impl_solve are called. More...
 
size_t n_impl_solve
 
size_t n_newton_iter
 
fstream output_file
 Output file. More...
 
double error
 error: For \( \nu=0 \), vdP reduces to linear oscillator, otherwise no analytic solution is available More...
 

Additional Inherited Members

- Protected Member Functions inherited from pfasst::encap::ImplicitSweeper< time >
void set_end_state ()
 
- Protected Attributes inherited from pfasst::encap::ImplicitSweeper< time >
Matrix< time > q_tilde
 
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_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::vdp::VdpSweeper< time >

Sweeper for the van der Pol oscillator.

\[ x' = y\\ y' = \nu*(1 - x^2)*y - x \]

Note that an analytical solution is available only for \( \nu=0 \), where the vdP simplifies to the standard linear oscillator. Hence, actual errors are computed only for \( \nu=0 \).

Definition at line 53 of file vdp_sweeper.hpp.

Member Typedef Documentation

template<typename time = pfasst::time_precision>
typedef encap::Encapsulation<time> pfasst::examples::vdp::VdpSweeper< time >::encap_type
private

Definition at line 57 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
typedef encap::VectorEncapsulation<double> pfasst::examples::vdp::VdpSweeper< time >::real_vector_type
private

Define a type for a complex PFASST vector encapsulation.

Definition at line 60 of file vdp_sweeper.hpp.

Constructor & Destructor Documentation

template<typename time = pfasst::time_precision>
pfasst::examples::vdp::VdpSweeper< time >::VdpSweeper ( double  nu,
double  x0,
double  y0 
)
inline

generic constructor; initialize all function call counters with zero.

Also open file to write solution.

Parameters
[in]coefficientnu of van-der-Pol oscillator

Definition at line 91 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
virtual pfasst::examples::vdp::VdpSweeper< time >::~VdpSweeper ( )
inlinevirtual

upon destruction, report final error and number of function calls and close output file.

Definition at line 110 of file vdp_sweeper.hpp.

References ML_LOG.

Member Function Documentation

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

computes and prints out error, which is meaningful only for \( \nu=0 \).

Also writes the current solution into an ascii file

Parameters
[in]tTime

Definition at line 124 of file vdp_sweeper.hpp.

References pfasst::examples::boris::max(), and ML_LOG.

+ Here is the call graph for this function:

template<typename time = pfasst::time_precision>
void pfasst::examples::vdp::VdpSweeper< time >::exact ( real_vector_type q,
time  t 
)
inline

if \( \nu=0 \), this computes the exact solution at time t.

For other values of \( \nu \), it just returns the inital value, because there is no analytical solution to evaluate.

Parameters
[in]tTime

For \( \nu=0 \) and given initial value \( x0 \), \( y0 \), the analytic solution reads

\[ x(t) = y0*sin(t) + x0*cos(t)\\ y(t) = -x0*sin(t) + y0*cos(t) \]

Note that for \( t=0 \), we recover \( x(0) = x0 \), \( y(0) = y0 \).

Definition at line 178 of file vdp_sweeper.hpp.

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

Definition at line 205 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
void pfasst::examples::vdp::VdpSweeper< 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.

Definition at line 214 of file vdp_sweeper.hpp.

References UNUSED.

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::get_errors ( )
inline

returns error, but does not update it!

Definition at line 144 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
void pfasst::examples::vdp::VdpSweeper< 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( u - \Delta t f(u) \right) u = b \) for \( u \) and set f_encap to \( f(u) \).

Solves the nonlinear equation

\[ P(q) := q - dt*f(q) - rhs = 0 \]

using Newton's method. For the vdP oscillator, it is

\[ P(q) = \left( x - dt*y - rhs_0 ; y - dt*\left( \nu*(1-x^2)*y -x \right) - rhs_1 \right) \]

The Newton update reads

\[ q_new = q + inv(P'(q))*(-P(q)) \]

Note that \( P'(q) = Id - dt*f'(q) \).

The inverse of \( P' \) has been computed symbolically, so that \( P'(q)*x \) can be evaluated directly.

The Jacobian of the right hand side of the van der Pol oscillator for \( q=[x;y] \) is

\[ Df(q) = \left( 0 1 ; -2*\nu*x*y - 1 \nu*(1-x^2) \right) \]

Because here we need to invert \( Id - dt*f(q) = b \), the corresponding Jacobian is

\[ J(q) = \left( 1 -dt ; dt*2*\nu*x*y+1 1 + dt*\nu(x^2-1) \right) \]

which can be inverted analytically

\[ inv(J(q)) = (1/c)*\left( dt*x^2-dt+1 dt ; -2*dt*\nu*x*y-dt 1\right) =: [a dt; b 1.0] \]

with \( c:=2*\nu*x*y*dt^2 + dt^2 + dt*x^2 - dt + 1 \)

Definition at line 234 of file vdp_sweeper.hpp.

References plot::dt, and UNUSED.

template<typename time = pfasst::time_precision>
void pfasst::examples::vdp::VdpSweeper< time >::post_predict ( )
inlineoverridevirtual

post prediction step, update error.

Reimplemented from pfasst::ISweeper< time >.

Definition at line 152 of file vdp_sweeper.hpp.

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

template<typename time = pfasst::time_precision>
void pfasst::examples::vdp::VdpSweeper< time >::post_sweep ( )
inlineoverridevirtual

post sweep, update error.

Reimplemented from pfasst::ISweeper< time >.

Definition at line 162 of file vdp_sweeper.hpp.

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

Member Data Documentation

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::error
private

error: For \( \nu=0 \), vdP reduces to linear oscillator, otherwise no analytic solution is available

Definition at line 81 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::vdp::VdpSweeper< time >::n_f_impl_eval
private

Counters for how often f_impl_eval and impl_solve are called.

Definition at line 75 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::vdp::VdpSweeper< time >::n_impl_solve
private

Definition at line 75 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::vdp::VdpSweeper< time >::n_newton_iter
private

Definition at line 75 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
size_t pfasst::examples::vdp::VdpSweeper< time >::newton_maxit
private

Maximum number of iterations for the nonlinear Newton solver.

Definition at line 69 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::newton_tol
private

Tolerance for the nonlinear Newton iteration.

Definition at line 72 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::nu
private

Parameter \( \nu \) in the van-der-Pol oscillator.

Definition at line 63 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
fstream pfasst::examples::vdp::VdpSweeper< time >::output_file
private

Output file.

Definition at line 78 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::x0
private

Starting values.

Definition at line 66 of file vdp_sweeper.hpp.

template<typename time = pfasst::time_precision>
double pfasst::examples::vdp::VdpSweeper< time >::y0
private

Definition at line 66 of file vdp_sweeper.hpp.


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