KASKADE 7 development version
|
#include <util.hh>
Public Types | |
using | field_type = typename Vector::value_type |
using | RealVector = Dune::DynamicVector< double > |
using | RealMatrix = Dune::DynamicMatrix< double > |
Public Member Functions | |
field_type | sdcContractionFactor (std::vector< Vector > const &yPrev, std::vector< Vector > const &yCurrent, std::vector< Vector > const &yNext, Norm &norm) |
virtual RealVector const & | computeAlphaVec (RealVector const &timePoints, RealMatrix const &integrationMatrix, std::vector< Vector > const &yPrev, std::vector< Vector > const &yCurrent, std::vector< Vector > const &yNext)=0 |
virtual | ~SDCUtil () |
using Kaskade::SDCUtil< Vector, Norm >::field_type = typename Vector::value_type |
using Kaskade::SDCUtil< Vector, Norm >::RealMatrix = Dune::DynamicMatrix<double> |
using Kaskade::SDCUtil< Vector, Norm >::RealVector = Dune::DynamicVector<double> |
|
inlinevirtual |
|
pure virtual |
The function definition changes with the associated norm. For derivation details for the vectors \( \alpha \) and \( \Gamma \) see M. Weiser, S. Ghosh: Adaptive inexact SDC Methods. The vector \( \alpha \) is given by \( \alpha = \{\alpha_1, \ldots , \alpha_N \}\), where
\begin{eqnarray*} \alpha_i & = & \sum_{k=1}^N |S_{ki}| + (1 - \delta_{i,N})\,(t_{i+1} - t_i)\,(1 + \rho), \end{eqnarray*}
this is computed in case of \( 1\)-norm and in the case of \( \max\)-norm we have \( \Gamma = \{\Gamma_1, \ldots, \Gamma_N \} \), where
\begin{eqnarray*} \Gamma_i & = & \max_{1\leq n \leq N} \alpha_{n,i}, \textrm{ such that } \\ \alpha_{n,i} & = & |S_{ni}| + \delta_{i,n-1}\,(t_i - t_{i-1})\,(1+\rho) \end{eqnarray*}
[in] | timePoints | : Time points generated using the points() method of SDCTimeGrid. |
[in] | integrationMatrix | : Integration matrix generated using the integrationMatrix() method of SDCTimeGrid. Computed once for a given time grid. |
[in] | yPrev | : Estimate of exact solution \( y \) using SDC iteration step |
[in] | yCurrent | : The following estimate of \( y \) using SDC iteration step |
[in] | yNext | : The next estimate of \( y \) after yCurrent using SDC iteration step |
Implemented in Kaskade::SDCUtilOneNorm< Vector, Norm >, and Kaskade::SDCUtilMaxNorm< Vector, Norm >.
Vector::value_type Kaskade::SDCUtil< Vector, Norm >::sdcContractionFactor | ( | std::vector< Vector > const & | yPrev, |
std::vector< Vector > const & | yCurrent, | ||
std::vector< Vector > const & | yNext, | ||
Norm & | norm | ||
) |
This function computes an estimate of the SDC contraction factor ( \( \rho \)) given three consecutive iterations of \( y \) and depends on the given norm. The estimate of \( \rho \) is given by:
\begin{eqnarray*} \rho & = & \frac{\|y^{[j+1]} - y^{[j]}\|}{\|y^{[j]} - y^{[j-1]}\| + \|y^{[j+1]} - y^{[j]}\|} \end{eqnarray*}
[in] | yPrev | a Vector denoting previous iterate of \( y \). |
[in] | yCurrent | a Vector denoting current iterate of \( y \). |
[in] | yNext | a Vector denoting next iterate of \( y \) as obtained from SDC iteration step. |
[in] | norm | an object of abstract base class type Norm, where we can provide the specific norm we consider for a problem. |
Definition at line 466 of file util.hh.
Referenced by Kaskade::SDCUtilOneNorm< Vector, Norm >::computeAlphaVec().