DUNE-ACFEM (unstable)

femschemeinterface.hh
Go to the documentation of this file.
1 
5 #ifndef __FEMSCHEME_INTERFACE_HH__
6 #define __FEMSCHEME_INTERFACE_HH__
7 
8 /*********************************************************/
9 
10 
11 #include <cstddef>
12 #include <string>
13 
14 namespace Dune {
15 
16  namespace ACFem {
17 
28  {
29  public:
31  virtual void initialize() = 0;
32 
34  virtual std::string name() const = 0;
35 
47  virtual void solve(bool forceMatrixAssembling = true) = 0;
48 
51  virtual int output() = 0;
52 
54  virtual bool converged() const = 0;
55 
57  virtual double residual() const = 0;
58 
60  virtual double error() const = 0;
61 
63  virtual size_t size() const = 0;
64  };
65 
69  : public BasicFemScheme
70  {
71  public:
73  virtual bool mark(const double tolerance) = 0;
74 
76  virtual double estimate() = 0;
77 
79  virtual void adapt() = 0;
80  };
81 
82  class TransientAdaptiveFemScheme
83  : public AdaptiveFemScheme
84  {
85  public:
87  virtual void restart(const std::string& name = "") = 0;
88 
90  virtual void next() = 0;
91 
93  virtual double timeEstimate() = 0;
94 
97  virtual double initialEstimate() = 0;
98 
100  virtual bool initialMarking(const double tolerance) = 0;
101 
102  };
103 
105 
106  } // ACFem::
107 
108 } // Dune::
109 
110 #endif // __FEMSCHEME_INTERFACE_HH__
Abstract space adaptative FEM scheme.
Definition: femschemeinterface.hh:70
virtual bool mark(const double tolerance)=0
mark elements for adaptation
virtual void adapt()=0
do the adaptation for a given marking
virtual double estimate()=0
calculate error estimator
Abstract non-adaptive basic FEM scheme.
Definition: femschemeinterface.hh:28
virtual size_t size() const =0
return some measure about the number of DOFs in use
virtual double error() const =0
calculate L2/H1 error
virtual int output()=0
data I/O, return -1 if no data has been written otherwise the sequence number of the file
virtual void solve(bool forceMatrixAssembling=true)=0
Solve the system.
virtual std::string name() const =0
name of the Fem Scheme
virtual void initialize()=0
initialize the solution
virtual bool converged() const =0
check whether solver has converged
virtual double residual() const =0
calculate residual (in small l^2)
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 2, 22:35, 2024)