LifeV
TimeAdvanceBDFVariableStep< FEVectorType > Class Template Reference

TimeAdvanceBDFVariableStep - Backward differencing formula time discretization with non-uniform time step. More...

#include <TimeAdvanceBDFVariableStep.hpp>

+ Collaboration diagram for TimeAdvanceBDFVariableStep< FEVectorType >:

Private Attributes

UInt M_order
 Order of the BDF derivative/extrapolation: the time-derivative coefficients vector has size n+1, the extrapolation vector has size n. More...
 
container_Type M_alpha
 Coefficients $ \alpha_i $ of the time bdf discretization. More...
 
container_Type M_beta
 Coefficients $ \beta_i $ of the extrapolation. More...
 
container_Type M_timeStep
 container_Type $ \delta_t $ of time intervals More...
 
container_Type M_timeStepBack
 
feVectorPtrContainer_Type M_unknowns
 Last n state vectors. More...
 
feVectorPtrContainer_Type M_unknownsBack
 

Public Types

typedef ScalarVector container_Type
 
typedef FEVectorType feVector_Type
 
typedef std::shared_ptr< feVector_TypefeVectorPtr_Type
 
typedef std::vector< feVectorPtr_TypefeVectorPtrContainer_Type
 

Constructor & Destructor

 TimeAdvanceBDFVariableStep ()
 Empty constructor. More...
 
 ~TimeAdvanceBDFVariableStep ()
 Destructor. More...
 

Methods

void setup (const UInt order)
 Setup the TimeAdvanceBDFVariableStep with order n. More...
 
void setInitialCondition (feVector_Type u0, Real const timeStep, bool startup=0)
 Initialize all the entries of the unknown vector to be derived with the vector u0 (duplicated if startup=0) More...
 
void setInitialCondition (std::vector< feVector_Type > u0s, Real const timeStep)
 Initialize all the entries of the unknown vector to be derived with a set of vectors u0s. More...
 
template<typename FunctionType , typename FESpaceType >
void setInitialCondition (const FunctionType &u0Function, feVector_Type &u0Vector, FESpaceType &feSpace, Real t0, Real timeStep)
 Initialize all the entries of the unknonwn vectors with a given function. More...
 
feVector_Type rhsContribution () const
 Returns the right hand side $ \bar{p} $ of the time derivative formula divided by dt. More...
 
feVector_Type rhsContributionTimeStep (Real timeStep=1) const
 Returns the right hand side $ \bar{p} $ of the time derivative formula divided by timeStep (backward compatibility version, will be discontinued) More...
 
feVector_Type extrapolateSolution () const
 Compute the polynomial extrapolation approximation of order n-1 of u^{n+1} defined by the n stored state vectors. More...
 
void shiftRight (feVector_Type const &uCurrent)
 Update the vectors of the previous time steps by shifting on the right the old values. More...
 
void shiftRight (feVector_Type const &uCurrent, Real timeStepNew)
 Update the vectors of the previous time steps by shifting on the right the old values. More...
 
void storeSolution ()
 Save the current vector M_unknowns and the current vector M_timeStep. More...
 
void restoreSolution ()
 Restore the vector M_unknowns and the vector M_timeStep with the ones saved with store_unk() More...
 
void restoreSolution (Real timeStep)
 It is equivalent to do : storeSolution() + setTimeStep() More...
 
void showMe () const
 Show informations about the BDF. More...
 

Set Methods

void setTimeStep (Real timeStep)
 Replace the current time-step with timeStep and computes the coefficients a_i and beta_i as functions of M_timeStep. More...
 

Get Methods

const RealcoefficientDerivative (UInt i) const
 Return the i-th coefficient of the time derivative alpha_i. More...
 
Real coefficientDerivativeOverTimeStep (UInt i) const
 Return the i-th coefficient of the time derivative alpha_i divided by dt. More...
 
const RealcoefficientExtrapolation (UInt i) const
 Return the i-th coefficient of the time extrapolation beta_i. More...
 
const container_TypetimeStepVector () const
 Return the vector of the time steps, ordered starting from the most recent one. More...
 
const feVectorPtrContainer_TypestateVector () const
 Return a vector with the last n state vectors. More...
 
const feVectorPtrContainer_TypestateVectorBack () const
 
void computeCoefficient ()
 Private methods. More...
 

Detailed Description

template<typename FEVectorType = VectorEpetra>
class LifeV::TimeAdvanceBDFVariableStep< FEVectorType >

TimeAdvanceBDFVariableStep - Backward differencing formula time discretization with non-uniform time step.

A differential equation of the form

$ M u' = A u + f $

is discretized in time as

$ M p'(t_{k+1}) = A u_{k+1} + f_{k+1} $

where p denotes the polynomial of order n in t that interpolates $ (t_i,u_i) $ for i = k-n+1,...,k+1.

The approximative time derivative $ p'(t_{k+1}) $ is a linear combination of state vectors $ u_i $:

$ p'(t_{k+1}) = \frac{1}{\Delta t} (\alpha_0 u_{k+1} - \sum_{i=0}^n \alpha_i u_{k+1-i} )$

Thus we have

$ \frac{\alpha_0}{\Delta t} M u_{k+1} = A u_{k+1} + f + M \bar{p} $

with

$ \bar{p} = \frac{1}{\Delta t} \sum_{i=1}^n \alpha_i u_{k+1-i} $

This class stores the n last state vectors in order to be able to calculate $ \bar{p} $. It also provides $ \alpha_i $ and can extrapolate the the new state from the n last states with a polynomial of order n-1:

$ u_{k+1} \approx \sum_{i=0}^{n-1} \beta_i u_{k-i} $

The class allows to change the time-steps, recomputing the coefficients. to change $ \Delta t $ at each time-step, use shift_right(uCurrent, timeStepNew) instead of shift_right(uCurrent).

Other methods that could be helpful while adapting in time are: set_deltat, store_unk and restore_unk

Definition at line 97 of file TimeAdvanceBDFVariableStep.hpp.

Member Typedef Documentation

◆ container_Type

Definition at line 103 of file TimeAdvanceBDFVariableStep.hpp.

◆ feVector_Type

typedef FEVectorType feVector_Type

Definition at line 104 of file TimeAdvanceBDFVariableStep.hpp.

◆ feVectorPtr_Type

typedef std::shared_ptr< feVector_Type > feVectorPtr_Type

Definition at line 105 of file TimeAdvanceBDFVariableStep.hpp.

◆ feVectorPtrContainer_Type

Definition at line 106 of file TimeAdvanceBDFVariableStep.hpp.

Constructor & Destructor Documentation

◆ TimeAdvanceBDFVariableStep()

Empty constructor.

Definition at line 332 of file TimeAdvanceBDFVariableStep.hpp.

◆ ~TimeAdvanceBDFVariableStep()

Destructor.

Definition at line 342 of file TimeAdvanceBDFVariableStep.hpp.

Member Function Documentation

◆ setup()

void setup ( const UInt  order)

Setup the TimeAdvanceBDFVariableStep with order n.

Parameters
orderorder of the BDF

Definition at line 349 of file TimeAdvanceBDFVariableStep.hpp.

◆ setInitialCondition() [1/3]

void setInitialCondition ( feVector_Type  u0,
Real const  timeStep,
bool  startup = 0 
)

Initialize all the entries of the unknown vector to be derived with the vector u0 (duplicated if startup=0)

If startup = 0, it initializes all the entries of the unknown vectors with a given function The array of initial conditions needed by the selected BDF is initialized as follows: unk=[ u0(t0), u0(t0-dt), u0(t0-2*dt), ...] When startUp = true, only M_unknown[0] is initialized (with u0). At the first step, class Bfd computes the coefficients of BDF1, at the second step, the coefficient of BDF2, and so on, until it reaches the wanted method. Second order of accuracy is obtained using this procedure to startup BDF2. Using the startup procedure with BDF3, the accuracy order is limited to the second order; anyway it is better to use the startup procedure instead of taking u{-2} = u_{-1} = u_0.

Parameters
u0initial condition vector
timeSteptime step
startup

Definition at line 370 of file TimeAdvanceBDFVariableStep.hpp.

◆ setInitialCondition() [2/3]

void setInitialCondition ( std::vector< feVector_Type u0s,
Real const  timeStep 
)

Initialize all the entries of the unknown vector to be derived with a set of vectors u0s.

Parameters
u0sinitial condition vectors
timeSteptime step

Definition at line 395 of file TimeAdvanceBDFVariableStep.hpp.

◆ setInitialCondition() [3/3]

void setInitialCondition ( const FunctionType &  u0Function,
feVector_Type u0Vector,
FESpaceType &  feSpace,
Real  t0,
Real  timeStep 
)

Initialize all the entries of the unknonwn vectors with a given function.

The array of initial conditions needed by the selected BDF is initialized as follows: M_unknown=[ u0Function(t0), u0Function(t0-dt), u0Function(t0-2*dt), ...] For the space dependence of the initial conditions we need informations on:

Parameters
u0Functionthe function we want to interpolate
u0Vectorcorrected mapped vector, it contains the u0Function(t0) in output
feSpacefinite element space
t0initial time (t0) and the
timeSteptime step (for solutions before the initial instant)

Based on the NavierStokesHandler::initialize by M. Fernandez

Definition at line 426 of file TimeAdvanceBDFVariableStep.hpp.

◆ rhsContribution()

TimeAdvanceBDFVariableStep< FEVectorType >::feVector_Type rhsContribution ( ) const

Returns the right hand side $ \bar{p} $ of the time derivative formula divided by dt.

Returns
a feVector_Type, by default feVector_Type = VectorEpetra

Definition at line 448 of file TimeAdvanceBDFVariableStep.hpp.

◆ rhsContributionTimeStep()

TimeAdvanceBDFVariableStep< FEVectorType >::feVector_Type rhsContributionTimeStep ( Real  timeStep = 1) const

Returns the right hand side $ \bar{p} $ of the time derivative formula divided by timeStep (backward compatibility version, will be discontinued)

Parameters
timeSteptime step, 1 by default
Returns
a feVector_Type, by default feVector_Type = VectorEpetra

Definition at line 464 of file TimeAdvanceBDFVariableStep.hpp.

◆ extrapolateSolution()

TimeAdvanceBDFVariableStep< FEVectorType >::feVector_Type extrapolateSolution ( ) const

Compute the polynomial extrapolation approximation of order n-1 of u^{n+1} defined by the n stored state vectors.

Definition at line 480 of file TimeAdvanceBDFVariableStep.hpp.

◆ shiftRight() [1/2]

void shiftRight ( feVector_Type const &  uCurrent)

Update the vectors of the previous time steps by shifting on the right the old values.

Parameters
uCurrentcurrent (new) value of the state vector

Definition at line 496 of file TimeAdvanceBDFVariableStep.hpp.

◆ shiftRight() [2/2]

void shiftRight ( feVector_Type const &  uCurrent,
Real  timeStepNew 
)

Update the vectors of the previous time steps by shifting on the right the old values.

Set the the new time step and shift right the old time steps.

Parameters
uCurrentcurrent (new) value of the state vector
timeStepNewnew value of the time step

Definition at line 521 of file TimeAdvanceBDFVariableStep.hpp.

◆ storeSolution()

void storeSolution ( )

Save the current vector M_unknowns and the current vector M_timeStep.

Definition at line 548 of file TimeAdvanceBDFVariableStep.hpp.

◆ restoreSolution() [1/2]

void restoreSolution ( )

Restore the vector M_unknowns and the vector M_timeStep with the ones saved with store_unk()

Definition at line 563 of file TimeAdvanceBDFVariableStep.hpp.

◆ restoreSolution() [2/2]

void restoreSolution ( Real  timeStep)

It is equivalent to do : storeSolution() + setTimeStep()

Parameters
timeStepnew time step

Definition at line 580 of file TimeAdvanceBDFVariableStep.hpp.

◆ showMe()

void showMe ( ) const

Show informations about the BDF.

Definition at line 598 of file TimeAdvanceBDFVariableStep.hpp.

◆ setTimeStep()

void setTimeStep ( Real  timeStep)

Replace the current time-step with timeStep and computes the coefficients a_i and beta_i as functions of M_timeStep.

Parameters
timeSteptime step

Definition at line 617 of file TimeAdvanceBDFVariableStep.hpp.

◆ coefficientDerivative()

const Real & coefficientDerivative ( UInt  i) const

Return the i-th coefficient of the time derivative alpha_i.

Parameters
iindex of the coefficient

Definition at line 628 of file TimeAdvanceBDFVariableStep.hpp.

◆ coefficientDerivativeOverTimeStep()

Real coefficientDerivativeOverTimeStep ( UInt  i) const

Return the i-th coefficient of the time derivative alpha_i divided by dt.

Parameters
iindex of the coefficient

Definition at line 639 of file TimeAdvanceBDFVariableStep.hpp.

◆ coefficientExtrapolation()

const Real & coefficientExtrapolation ( UInt  i) const

Return the i-th coefficient of the time extrapolation beta_i.

Parameters
iindex of the coefficient

Definition at line 650 of file TimeAdvanceBDFVariableStep.hpp.

◆ timeStepVector()

const container_Type& timeStepVector ( ) const
inline

Return the vector of the time steps, ordered starting from the most recent one.

Definition at line 261 of file TimeAdvanceBDFVariableStep.hpp.

◆ stateVector()

const feVectorPtrContainer_Type& stateVector ( ) const
inline

Return a vector with the last n state vectors.

Definition at line 267 of file TimeAdvanceBDFVariableStep.hpp.

◆ stateVectorBack()

const feVectorPtrContainer_Type& stateVectorBack ( ) const
inline

Definition at line 272 of file TimeAdvanceBDFVariableStep.hpp.

◆ computeCoefficient()

void computeCoefficient ( )
private

Private methods.

Computes the coefficients a_i and beta_i as functions of _M_delta_t.

Arbitrary order, variable time step BDF coefficients. Reference: Comincioli pag. 618-619 Matlab implementation: function [alpha0, alpha, beta] = compute_coef(dts, order) assert( (length(dts) == order) ); rho = dts(1)./cumsum(dts);

alpha0 = sum(rho); alpha = zeros(order,1); beta = zeros(order,1);

for j = 1:order ind = setdiff([1:order],j); beta(j) = 1/prod(1-rho(ind)/rho(j)); alpha(j) = rho(j)*beta(j); end

end

Definition at line 662 of file TimeAdvanceBDFVariableStep.hpp.

+ Here is the caller graph for this function:

Field Documentation

◆ M_order

UInt M_order
private

Order of the BDF derivative/extrapolation: the time-derivative coefficients vector has size n+1, the extrapolation vector has size n.

Definition at line 310 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_alpha

container_Type M_alpha
private

Coefficients $ \alpha_i $ of the time bdf discretization.

Definition at line 313 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_beta

container_Type M_beta
private

Coefficients $ \beta_i $ of the extrapolation.

Definition at line 316 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_timeStep

container_Type M_timeStep
private

container_Type $ \delta_t $ of time intervals

Definition at line 319 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_timeStepBack

container_Type M_timeStepBack
private

Definition at line 320 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_unknowns

feVectorPtrContainer_Type M_unknowns
private

Last n state vectors.

Definition at line 323 of file TimeAdvanceBDFVariableStep.hpp.

◆ M_unknownsBack

feVectorPtrContainer_Type M_unknownsBack
private

Definition at line 324 of file TimeAdvanceBDFVariableStep.hpp.


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