Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
Algorithms::ThetaTimestepping< VectorType > Class Template Reference

#include <deal.II/algorithms/theta_timestepping.h>

Inheritance diagram for Algorithms::ThetaTimestepping< VectorType >:
[legend]

Public Member Functions

 ThetaTimestepping (OperatorBase &op_explicit, OperatorBase &op_implicit)
 
virtual void operator() (AnyData &out, const AnyData &in) override
 
virtual void notify (const Event &) override
 
void set_output (OutputOperator< VectorType > &output)
 
void parse_parameters (ParameterHandler &param)
 
double current_time () const
 
double theta () const
 
double theta (double new_theta)
 
const TimestepDataexplicit_data () const
 
const TimestepDataimplicit_data () const
 
TimestepControltimestep_control ()
 
- Public Member Functions inherited from Algorithms::OperatorBase
virtual ~OperatorBase () override=default
 
void clear_events ()
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &param)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Attributes

TimestepControl control
 
double vtheta
 
bool adaptive
 
TimestepData d_explicit
 
TimestepData d_implicit
 
SmartPointer< OperatorBase, ThetaTimestepping< VectorType > > op_explicit
 
SmartPointer< OperatorBase, ThetaTimestepping< VectorType > > op_implicit
 
SmartPointer< OutputOperator< VectorType >, ThetaTimestepping< VectorType > > output
 

Additional Inherited Members

- Protected Attributes inherited from Algorithms::OperatorBase
Event notifications
 

Detailed Description

template<typename VectorType>
class Algorithms::ThetaTimestepping< VectorType >

Application class performing the theta timestepping scheme.

The theta scheme is an abstraction of implicit and explicit Euler schemes, the Crank-Nicholson scheme and linear combinations of those. The choice of the actual scheme is controlled by the parameter theta as follows.

For fixed theta, the Crank-Nicholson scheme is the only second order scheme. Nevertheless, further stability may be achieved by choosing theta larger than ½, thereby introducing a first order error term. In order to avoid a loss of convergence order, the adaptive theta scheme can be used, where theta=½+c dt.

Assume that we want to solve the equation u' + F(u) = 0 with a step size k. A step of the theta scheme can be written as

\[ M u_{n+1} + \theta k F(u_{n+1}) = M u_n - (1-\theta)k F(u_n). \]

Here, M is the mass matrix. We see, that the right hand side amounts to an explicit Euler step with modified step size in weak form (up to inversion of M). The left hand side corresponds to an implicit Euler step with modified step size (right hand side given). Thus, the implementation of the theta scheme will use two Operator objects, one for the explicit, one for the implicit part. Each of these will use its own TimestepData to account for the modified step sizes (and different times if the problem is not autonomous). Note that once the explicit part has been computed, the left hand side actually constitutes a linear or nonlinear system which has to be solved.

Usage AnyData

ThetaTimestepping uses AnyData for communicating vectors and time step information. With outer or inner Operator objects. It does not use itself the input vectors provided, but forwards them to the explicit and implicit operators.

Vector data

The explicit Operator op_explicit receives in its input in first place the vector "Previous iterate", which is the solution value after the previous timestep. It is followed by all vectors provided to ThetaTimestepping::operator() as input argument. op_explicit is supposed to write its result into the first position of its output argument, labeled "Result".

The implicit Operator op_implicit receives the result of op_explicit in its first input vector labeled "Previous time". It is followed by all vectors provided to ThetaTimestepping::operator() as input argument. The output of op_implicit is directly written into the output argument given to ThetaTimestepping.

Scalar data

Since the introduction of AnyData, ThetaTimestepping is able to communicate the current time step information through AnyData as well. Therefore, the AnyData objects handed as input to op_explicit and op_implicit contain two entries of type const double* named "Time" and "Timestep". Note that "Time" refers to the time at the beginning of the current step for op_explicit and at the end for op_implicit, respectively.

Usage of ThetaTimestepping

The use ThetaTimestepping is more complicated than for instance Newton, since the inner operators will usually need to access the TimeStepData. Thus, we have a circular dependency of information, and we include the following example for its use. It can be found in examples/doxygen/theta_timestepping.cc

First, we define the two operators used by ThetaTimestepping and call them Implicit and Explicit. They both share the public interface of Operator, and additionally provide storage for the matrices to be used and a pointer to TimestepData. Note that we do not use a SmartPointer here, since the TimestepData will be destroyed before the operator.

These operators will be implemented after the main program. But let us look first at how they get used. First, let us define a matrix to be used for our system and also an OutputOperator in order to write the data of each timestep to a file.

Now we create objects for the implicit and explicit parts of the steps as well as the ThetaTimestepping itself. We initialize the timestepping with the output operator in order to be able to see the output in every step.

The next step is providing the vectors to be used. value is filled with the initial value and is also the vector where the solution at each timestep will be. Because the interface of Operator has to be able to handle several vectors, we need to store it in an AnyData object. Since our problem has no additional parameters, the input AnyData object remains empty.

Finally, we are ready to tell the solver, that we are starting at the initial timestep and run it.

First the constructor, which simply copies the system matrix into the member pointer for later use.

Now we need to study the application of the implicit and explicit operator. We assume that the pointer matrix points to the matrix created in the main program (the constructor did this for us). Here, we first get the time step size from the AnyData object that was provided as input. Then, if we are in the first step or if the timestep has changed, we fill the local matrix \(m\), such that with the given matrix \(M\), it becomes

\[ m = I - \Delta t M. \]

After we have worked off the notifications, we clear them, such that the matrix is only generated when necessary.

Now we multiply the input vector with the new matrix and store on output.

that we change the sign in front of the timestep and use the inverse of t he matrix.

Author
Guido Kanschat
Date
2010

Definition at line 191 of file theta_timestepping.h.

Constructor & Destructor Documentation

◆ ThetaTimestepping()

template<typename VectorType >
Algorithms::ThetaTimestepping< VectorType >::ThetaTimestepping ( OperatorBase op_explicit,
OperatorBase op_implicit 
)

Constructor, receiving the two operators stored in op_explicit and op_implicit. For their meaning, see the description of those variables.

Member Function Documentation

◆ operator()()

template<typename VectorType >
virtual void Algorithms::ThetaTimestepping< VectorType >::operator() ( AnyData out,
const AnyData in 
)
overridevirtual

The timestepping scheme.

Parameters
inis ignored by ThetaTimestepping, but is merged into the AnyData objects used as input for the operators op_explicit and op_implicit.
outin its first argument must contain a pointer to a VectorType instance, which contains the initial value when the operator is called. It contains the final value when the operator returns.

Implements Algorithms::OperatorBase.

◆ notify()

template<typename VectorType >
virtual void Algorithms::ThetaTimestepping< VectorType >::notify ( const Event )
overridevirtual

Register an event triggered by an outer iteration.

Reimplemented from Algorithms::OperatorBase.

◆ set_output()

template<typename VectorType >
void Algorithms::ThetaTimestepping< VectorType >::set_output ( OutputOperator< VectorType > &  output)
inline

Define an operator which will output the result in each step. Note that no output will be generated without this.

Definition at line 367 of file theta_timestepping.h.

◆ declare_parameters()

template<typename VectorType >
static void Algorithms::ThetaTimestepping< VectorType >::declare_parameters ( ParameterHandler param)
static

Declare parameters in a parameter handler.

◆ parse_parameters()

template<typename VectorType >
void Algorithms::ThetaTimestepping< VectorType >::parse_parameters ( ParameterHandler param)

Read the parameters in the ParameterHandler.

◆ current_time()

template<typename VectorType >
double Algorithms::ThetaTimestepping< VectorType >::current_time ( ) const
inline

The current time in the timestepping scheme.

Definition at line 393 of file theta_timestepping.h.

◆ theta() [1/2]

template<typename VectorType >
double Algorithms::ThetaTimestepping< VectorType >::theta ( ) const
inline

The weight between implicit and explicit part.

Definition at line 375 of file theta_timestepping.h.

◆ theta() [2/2]

template<typename VectorType >
double Algorithms::ThetaTimestepping< VectorType >::theta ( double  new_theta)
inline

Set a new weight and return the old

Definition at line 383 of file theta_timestepping.h.

◆ explicit_data()

template<typename VectorType >
const TimestepData & Algorithms::ThetaTimestepping< VectorType >::explicit_data ( ) const
inline

The data handed to the op_explicit time stepping operator.

The time in here is the time at the beginning of the current step, the time step is (1-theta) times the actual time step.

Definition at line 344 of file theta_timestepping.h.

◆ implicit_data()

template<typename VectorType >
const TimestepData & Algorithms::ThetaTimestepping< VectorType >::implicit_data ( ) const
inline

The data handed to the op_implicit time stepping operator.

The time in here is the time at the beginning of the current step, the time step is theta times the actual time step.

Definition at line 352 of file theta_timestepping.h.

◆ timestep_control()

template<typename VectorType >
TimestepControl & Algorithms::ThetaTimestepping< VectorType >::timestep_control ( )
inline

Allow access to the control object.

Definition at line 360 of file theta_timestepping.h.

Member Data Documentation

◆ control

template<typename VectorType >
TimestepControl Algorithms::ThetaTimestepping< VectorType >::control
private

The object controlling the time step size and computing the new time in each step.

Definition at line 287 of file theta_timestepping.h.

◆ vtheta

template<typename VectorType >
double Algorithms::ThetaTimestepping< VectorType >::vtheta
private

The control parameter theta in the range [0,1]. It defaults to 0.5.

Definition at line 293 of file theta_timestepping.h.

◆ adaptive

template<typename VectorType >
bool Algorithms::ThetaTimestepping< VectorType >::adaptive
private

Use adaptive theta if true. Not yet implemented.

Definition at line 297 of file theta_timestepping.h.

◆ d_explicit

template<typename VectorType >
TimestepData Algorithms::ThetaTimestepping< VectorType >::d_explicit
private

The data for the explicit part of the scheme.

Definition at line 302 of file theta_timestepping.h.

◆ d_implicit

template<typename VectorType >
TimestepData Algorithms::ThetaTimestepping< VectorType >::d_implicit
private

The data for the implicit part of the scheme.

Definition at line 307 of file theta_timestepping.h.

◆ op_explicit

template<typename VectorType >
SmartPointer<OperatorBase, ThetaTimestepping<VectorType> > Algorithms::ThetaTimestepping< VectorType >::op_explicit
private

The operator computing the explicit part of the scheme. This will receive in its input data the value at the current time with name "Current time solution". It should obtain the current time and time step size from explicit_data().

Its return value is \( Mu+cF(u) \), where \(u\) is the current state vector, \(M\) the mass matrix, \(F\) the operator in space and \(c\) is the adjusted time step size \((1-\theta) \Delta t\).

Definition at line 320 of file theta_timestepping.h.

◆ op_implicit

template<typename VectorType >
SmartPointer<OperatorBase, ThetaTimestepping<VectorType> > Algorithms::ThetaTimestepping< VectorType >::op_implicit
private

The operator solving the implicit part of the scheme. It will receive in its input data the vector "Previous time". Information on the timestep should be obtained from implicit_data().

Its return value is the solution u of Mu-cF(u)=f, where f is the dual space vector found in the "Previous time" entry of the input data, M the mass matrix, F the operator in space and c is the adjusted time step size \( \theta \Delta t\)

Definition at line 332 of file theta_timestepping.h.

◆ output

template<typename VectorType >
SmartPointer<OutputOperator<VectorType>, ThetaTimestepping<VectorType> > Algorithms::ThetaTimestepping< VectorType >::output
private

The operator writing the output in each time step

Definition at line 338 of file theta_timestepping.h.


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