1#ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICIMPLICIT_HH 
    2#define DUNE_FEM_SOLVER_RUNGEKUTTA_BASICIMPLICIT_HH 
   16#include <dune/fem/solver/odesolverinterface.hh> 
   24  struct NoImplicitRungeKuttaSourceTerm
 
   27    bool operator() ( 
double time, 
double timeStepSize, 
int stage, 
const T &u, 
const std::vector< T * > &update, T &source )
 
   33    void limit( T& update, 
const double time ) {}
 
   36    double initialTimeStepEstimate ( 
double time, 
const T &u )
 const 
   42    double timeStepEstimate ()
 const 
   51  template< 
class HelmholtzOperator, 
class NonlinearSolver, 
class TimeStepControl, 
class SourceTerm = NoImplicitRungeKuttaSourceTerm >
 
   62    typedef HelmholtzOperator HelmholtzOperatorType;
 
   63    typedef NonlinearSolver NonlinearSolverType;
 
   64    typedef TimeStepControl TimeStepControlType;
 
   65    typedef SourceTerm SourceTermType;
 
   69    typedef typename TimeStepControlType::ParameterType ParameterType;
 
   70    typedef typename NonlinearSolver::ParameterType     NonlinearSolverParameterType;
 
   79    template< 
class ButcherTable >
 
   81                                    const ButcherTable &butcherTable,
 
   82                                    const TimeStepControlType &timeStepControl,
 
   83                                    const SourceTermType &sourceTerm,
 
   84                                    const NonlinearSolverParameterType& parameters )
 
   85    : helmholtzOp_( helmholtzOp ),
 
   86      nonlinearSolver_( parameters ),
 
   87      timeStepControl_( timeStepControl ),
 
   88      sourceTerm_( sourceTerm ),
 
   89      stages_( butcherTable.
stages() ),
 
   90      alpha_( butcherTable.A() ),
 
   93      c_( butcherTable.c() ),
 
   94      rhs_( 
"RK rhs", helmholtzOp_.space() ),
 
   96      update_( 
stages(), nullptr )
 
   98      setup( butcherTable );
 
  107    template< 
class ButcherTable >
 
  109                                    const ButcherTable &butcherTable,
 
  110                                    const TimeStepControlType &timeStepControl,
 
  111                                    const NonlinearSolverParameterType& parameters )
 
  112    : helmholtzOp_( helmholtzOp ),
 
  113      nonlinearSolver_( parameters ),
 
  114      timeStepControl_( timeStepControl ),
 
  115      stages_( butcherTable.
stages() ),
 
  116      alpha_( butcherTable.A() ),
 
  119      c_( butcherTable.c() ),
 
  120      rhs_( 
"RK rhs", helmholtzOp_.space() ),
 
  122      update_( 
stages(), nullptr )
 
  124      setup( butcherTable );
 
  127    template< 
class ButcherTable >
 
  129                                    const ButcherTable &butcherTable,
 
  130                                    const TimeStepControlType &timeStepControl,
 
  131                                    const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
 
  133        NonlinearSolverParameterType( parameter ) )
 
  137    template< 
class ButcherTable >
 
  139                                    const ButcherTable &butcherTable,
 
  140                                    const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
 
  142        TimeStepControlType( parameter ), NonlinearSolverParameterType( parameter ) )
 
  146    template< 
class ButcherTable >
 
  147    void setup( 
const ButcherTable& butcherTable )
 
  150      update_.resize( 
stages(), 
nullptr );
 
  151      updateStorage_.resize( 
stages() );
 
  154      for( 
int i = 0; i < 
stages(); ++i )
 
  156        std::ostringstream name;
 
  157        name << 
"RK stage " << i;
 
  158        updateStorage_[ i ].reset( 
new DestinationType( name.str(), helmholtzOp_.space() ) );
 
  159        update_[ i ] = updateStorage_[ i ].operator ->();
 
  164      for( 
int i = 0; i < 
stages(); ++i )
 
  166        gamma_[ i ] = AL[ i ][ i ];
 
  171      alpha_.
mtv( butcherTable.b(), beta_ );
 
  174      for( 
int i = 0; i < 
stages(); ++i )
 
  175        alpha_[ i ][ i ] = gamma_[ i ];
 
  177      for( 
int i = 0; i < 
stages(); ++i )
 
  180        for( 
int j = 0; j < i; ++j )
 
  181          gamma_[ i ] -= alpha_[ i ][ j ];
 
  185      for( 
int i = 0; i < 
stages(); ++i )
 
  186        delta_ -= beta_[ i ];
 
  193      const double time = timeStepControl_.time();
 
  195      helmholtzOp_.setTime( time );
 
  196      helmholtzOp_.initializeTimeStepSize( U0 );
 
  197      const double helmholtzEstimate = helmholtzOp_.timeStepEstimate();
 
  199      double sourceTermEstimate = sourceTerm_.initialTimeStepEstimate( time, U0 );
 
  201      if( sourceTermEstimate < 0.0 ) sourceTermEstimate = helmholtzEstimate ;
 
  203      timeStepControl_.initialTimeStepSize( helmholtzEstimate, sourceTermEstimate );
 
  209    void solve ( DestinationType &U, MonitorType &monitor )
 
  213      const double time = timeStepControl_.time();
 
  214      const double timeStepSize = timeStepControl_.timeStepSize();
 
  215      assert( timeStepSize > 0.0 );
 
  216      for( 
int s = 0; s < 
stages(); ++s )
 
  218        assert( update_[ s ] );
 
  220        DestinationType& updateStage = *update_[ s ];
 
  223        updateStage.assign( U );
 
  224        updateStage *= gamma_[ s ];
 
  225        for( 
int k = 0; k < s; ++k )
 
  226          updateStage.
axpy( alpha_[ s ][ k ], *update_[ k ] );
 
  228        const double stageTime = time + c_[ s ]*timeStepSize;
 
  229        if( sourceTerm_( time, timeStepSize, s, U, update_, rhs_ ) )
 
  231          updateStage.
axpy( alpha_[ s ][ s ]*timeStepSize, rhs_ );
 
  232          sourceTerm_.limit( updateStage, stageTime );
 
  236        helmholtzOp_.setTime( stageTime );
 
  237        helmholtzOp_.setLambda( 0 );
 
  238        helmholtzOp_( updateStage, rhs_ );
 
  241        helmholtzOp_.setLambda( alpha_[ s ][ s ]*timeStepSize );
 
  242        nonlinearSolver_.bind( helmholtzOp_ );
 
  243        nonlinearSolver_( rhs_, updateStage );
 
  244        nonlinearSolver_.unbind();
 
  247        monitor.newtonIterations_       += nonlinearSolver_.iterations();
 
  248        monitor.linearSolverIterations_ += nonlinearSolver_.linearIterations();
 
  251        if( !nonlinearSolver_.converged() )
 
  252          return timeStepControl_.reduceTimeStep( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
 
  256      if( timeStepControl_.computeError() )
 
  259        DestinationType Uerr( U );
 
  263        for( 
int s = 0; s < 
stages(); ++s )
 
  264          U.axpy( beta_[ s ], *update_[ s ] );
 
  267        Uerr.axpy( -1.0, U );
 
  268        const double errorU = Uerr.scalarProductDofs( Uerr );
 
  269        const double normU = U.scalarProductDofs( U );
 
  271        if( normU > 0 && errorU > 0 )
 
  273          error = std::sqrt( errorU / normU );
 
  275        std::cout << std::scientific << 
"Error in RK = " << error << 
" norm " << errorU << 
" " << normU << std::endl;
 
  282        for( 
int s = 0; s < 
stages(); ++s )
 
  283          U.axpy( beta_[ s ], *update_[ s ] );
 
  286      monitor.error_ = error;
 
  289      timeStepControl_.timeStepEstimate( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
 
  296      out << 
"Generic " << 
stages() << 
"-stage implicit Runge-Kutta solver.\\\\" << std::endl;
 
  300    double infNorm(
const DestinationType& U, 
const DestinationType& Uerr )
 const 
  302      typedef typename DestinationType :: ConstDofIteratorType ConstDofIteratorType ;
 
  303      const ConstDofIteratorType uend = U.dend();
 
  305      for( ConstDofIteratorType u = U.dbegin(), uerr = Uerr.dbegin(); u != uend; ++u, ++uerr )
 
  308        double uerrval = *uerr ;
 
  309        double div = std::abs( std::max( uval, uerrval ) );
 
  311        double norm = std::abs( uval - uerrval );
 
  312        if( std::abs(div) > 1e-12 )
 
  314        res = std::max( res, norm );
 
  319    HelmholtzOperatorType&   helmholtzOp_;
 
  320    NonlinearSolverType      nonlinearSolver_;
 
  321    TimeStepControl          timeStepControl_;
 
  322    SourceTerm               sourceTerm_;
 
  329    DestinationType rhs_;
 
  330    std::vector< std::unique_ptr< DestinationType > > updateStorage_;
 
  331    std::vector< DestinationType* > update_;
 
Implicit RungeKutta ODE solver.
Definition: basicimplicit.hh:54
 
void description(std::ostream &out) const
print description of ODE solver to out stream
Definition: basicimplicit.hh:294
 
int stages() const
return stages of RK solver (-1 if not implemented)
Definition: basicimplicit.hh:292
 
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const SourceTermType &sourceTerm, const NonlinearSolverParameterType ¶meters)
constructor
Definition: basicimplicit.hh:80
 
void initialize(const DestinationType &U0)
apply operator once to get dt estimate
Definition: basicimplicit.hh:191
 
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const NonlinearSolverParameterType ¶meters)
constructor
Definition: basicimplicit.hh:108
 
void solve(DestinationType &U, MonitorType &monitor)
solve the system
Definition: basicimplicit.hh:209
 
Interface class for ODE Solver.
Definition: odesolverinterface.hh:21
 
Monitor MonitorType
monitor type
Definition: odesolverinterface.hh:59
 
virtual void solve(DestinationType &u)
solve  where  is the internal operator.
Definition: odesolverinterface.hh:75
 
DestinationImp DestinationType
type of destination
Definition: odesolverinterface.hh:62
 
void invert(bool doPivoting=true)
Compute inverse.
 
MAT & leftmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the left to this matrix.
Definition: densematrix.hh:632
 
constexpr void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:392
 
constexpr derived_type & axpy(const field_type &a, const DenseVector< Other > &x)
vector space axpy operation ( *this += a x )
Definition: densevector.hh:576
 
general base for time providers
Definition: timeprovider.hh:36
 
This file implements a dense matrix with dynamic numbers of rows and columns.
 
This file implements a dense vector with a dynamic size.
 
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485