dune-pdelab  2.4-dev
default/localassembler.hh
Go to the documentation of this file.
1 #ifndef DUNE_PDELAB_DEFAULT_LOCAL_ASSEMBLER_HH
2 #define DUNE_PDELAB_DEFAULT_LOCAL_ASSEMBLER_HH
3 
4 #include <dune/typetree/typetree.hh>
5 
12 
13 namespace Dune{
14  namespace PDELab{
15 
30  template<typename GO, typename LOP, bool nonoverlapping_mode = false>
32  public Dune::PDELab::LocalAssemblerBase<typename GO::Traits::MatrixBackend,
33  typename GO::Traits::TrialGridFunctionSpaceConstraints,
34  typename GO::Traits::TestGridFunctionSpaceConstraints>
35  {
36 
37  // The GridOperator has to be a friend to modify the do{Pre,Post}Processing flags
38  template<typename, typename, typename,
39  typename, typename, typename, typename,
40  typename, typename, bool>
41  friend class GridOperator;
42 
43  public:
44 
47 
49  typedef typename Traits::Residual::ElementType RangeField;
50  typedef RangeField Real;
51 
54 
57 
60 
62  typedef typename GFSU::Traits::GridViewType GridView;
63 
65  typedef LOP LocalOperator;
66 
67  static const bool isNonOverlapping = nonoverlapping_mode;
68 
71  // Types of local function spaces
76 
78 
85 
91 
93  DefaultLocalAssembler (LOP & lop_, shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
94  : lop(lop_), weight(1.0), doPreProcessing(true), doPostProcessing(true),
95  pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this), jacobian_apply_engine(*this)
96  , _reconstruct_border_entries(isNonOverlapping)
97  {}
98 
100  DefaultLocalAssembler (LOP & lop_, const CU& cu_, const CV& cv_,
101  shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
102  : Base(cu_, cv_),
103  lop(lop_), weight(1.0), doPreProcessing(true), doPostProcessing(true),
104  pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this), jacobian_apply_engine(*this)
105  , _reconstruct_border_entries(isNonOverlapping)
106  {}
107 
111  void setTime(Real time_){
112  lop.setTime(time_);
113  }
114 
116  void setWeight(RangeField weight_){
117  weight = weight_;
118  }
119 
122  void preStage (Real time_, int r_) { lop.preStage(time_,r_); }
123  void preStep (Real time_, Real dt_, std::size_t stages_){ lop.preStep(time_,dt_,stages_); }
124  void postStep (){ lop.postStep(); }
125  void postStage (){ lop.postStage(); }
126  Real suggestTimestep (Real dt) const{return lop.suggestTimestep(dt); }
128 
130  {
131  return _reconstruct_border_entries;
132  }
133 
136 
139  LocalPatternAssemblerEngine & localPatternAssemblerEngine
141  {
142  pattern_engine.setPattern(p);
143  return pattern_engine;
144  }
145 
148  LocalResidualAssemblerEngine & localResidualAssemblerEngine
149  (typename Traits::Residual & r, const typename Traits::Solution & x)
150  {
151  residual_engine.setResidual(r);
152  residual_engine.setSolution(x);
153  return residual_engine;
154  }
155 
158  LocalJacobianAssemblerEngine & localJacobianAssemblerEngine
159  (typename Traits::Jacobian & a, const typename Traits::Solution & x)
160  {
161  jacobian_engine.setJacobian(a);
162  jacobian_engine.setSolution(x);
163  return jacobian_engine;
164  }
165 
168  LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine
169  (typename Traits::Residual & r, const typename Traits::Solution & x)
170  {
171  jacobian_apply_engine.setResidual(r);
172  jacobian_apply_engine.setSolution(x);
173  return jacobian_apply_engine;
174  }
175 
177 
182  static bool doAlphaVolume() { return LOP::doAlphaVolume; }
183  static bool doLambdaVolume() { return LOP::doLambdaVolume; }
184  static bool doAlphaSkeleton() { return LOP::doAlphaSkeleton; }
185  static bool doLambdaSkeleton() { return LOP::doLambdaSkeleton; }
186  static bool doAlphaBoundary() { return LOP::doAlphaBoundary; }
187  static bool doLambdaBoundary() { return LOP::doLambdaBoundary; }
188  static bool doAlphaVolumePostSkeleton() { return LOP::doAlphaVolumePostSkeleton; }
189  static bool doLambdaVolumePostSkeleton() { return LOP::doLambdaVolumePostSkeleton; }
190  static bool doSkeletonTwoSided() { return LOP::doSkeletonTwoSided; }
191  static bool doPatternVolume() { return LOP::doPatternVolume; }
192  static bool doPatternSkeleton() { return LOP::doPatternSkeleton; }
193  static bool doPatternBoundary() { return LOP::doPatternBoundary; }
194  static bool doPatternVolumePostSkeleton() { return LOP::doPatternVolumePostSkeleton; }
196 
201  void preProcessing(bool v)
202  {
203  doPreProcessing = v;
204  }
205 
210  void postProcessing(bool v)
211  {
212  doPostProcessing = v;
213  }
214 
215  private:
216 
218  LOP & lop;
219 
221  RangeField weight;
222 
225  bool doPreProcessing;
226 
229  bool doPostProcessing;
230 
233  LocalPatternAssemblerEngine pattern_engine;
234  LocalResidualAssemblerEngine residual_engine;
235  LocalJacobianAssemblerEngine jacobian_engine;
236  LocalJacobianApplyAssemblerEngine jacobian_apply_engine;
238 
239  bool _reconstruct_border_entries;
240 
241  };
242 
243  }
244 }
245 #endif
GO::Traits::Range Residual
The type of the range (residual).
Definition: assemblerutilities.hh:54
static bool doAlphaSkeleton()
Definition: default/localassembler.hh:184
static bool doAlphaBoundary()
Definition: default/localassembler.hh:186
static bool doLambdaVolume()
Definition: default/localassembler.hh:183
void preStage(Real time_, int r_)
Definition: default/localassembler.hh:122
Create a local function space from a global function space.
Definition: localfunctionspace.hh:668
GO::Traits::Jacobian Jacobian
The type of the jacobian.
Definition: assemblerutilities.hh:61
static bool doSkeletonTwoSided()
Definition: default/localassembler.hh:190
Dune::PDELab::LocalFunctionSpace< GFSU, Dune::PDELab::TrialSpaceTag > LFSU
Definition: default/localassembler.hh:72
GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace
The trial grid function space.
Definition: assemblerutilities.hh:26
MatrixBackend::template Pattern< Jacobian, TestGridFunctionSpace, TrialGridFunctionSpace > MatrixPattern
The matrix pattern.
Definition: assemblerutilities.hh:68
Traits::TrialGridFunctionSpace GFSU
Definition: default/localassembler.hh:52
void preProcessing(bool v)
Definition: default/localassembler.hh:201
The local assembler for DUNE grids.
Definition: default/localassembler.hh:31
void setPattern(Pattern &pattern_)
Definition: default/patternengine.hh:92
void setJacobian(Jacobian &jacobian_)
Definition: default/jacobianengine.hh:106
static bool doPatternVolume()
Definition: default/localassembler.hh:191
DefaultLocalResidualAssemblerEngine< DefaultLocalAssembler > LocalResidualAssemblerEngine
Definition: default/localassembler.hh:82
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: default/localassembler.hh:159
GFSU::Traits::GridViewType GridView
The current grid view type.
Definition: default/localassembler.hh:62
DefaultLocalAssembler(LOP &lop_, shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
Constructor with empty constraints.
Definition: default/localassembler.hh:93
static bool doPatternVolumePostSkeleton()
Definition: default/localassembler.hh:194
Standard grid operator implementation.
Definition: gridoperator.hh:34
void setTime(Real time_)
Definition: default/localassembler.hh:111
void setSolution(const Solution &solution_)
Definition: default/jacobianengine.hh:115
static bool doLambdaVolumePostSkeleton()
Definition: default/localassembler.hh:189
GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints
The type of the trial grid function space constraints.
Definition: assemblerutilities.hh:33
LFSIndexCache< LFSV, CV > LFSVCache
Definition: default/localassembler.hh:75
Traits::TestGridFunctionSpaceConstraints CV
Definition: default/localassembler.hh:56
Traits::TrialGridFunctionSpaceConstraints CU
Definition: default/localassembler.hh:55
LFSIndexCache< LFSU, CU > LFSUCache
Definition: default/localassembler.hh:74
Dune::PDELab::LocalAssemblerTraits< GO > Traits
The traits class.
Definition: default/localassembler.hh:46
void postStage()
Definition: default/localassembler.hh:125
bool reconstructBorderEntries() const
Definition: default/localassembler.hh:129
LOP LocalOperator
The local operator.
Definition: default/localassembler.hh:65
static bool doPatternBoundary()
Definition: default/localassembler.hh:193
void setWeight(RangeField weight_)
Notifies the assembler about the current weight of assembling.
Definition: default/localassembler.hh:116
Base class for local assembler.
Definition: assemblerutilities.hh:210
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: default/localassembler.hh:149
Real suggestTimestep(Real dt) const
Definition: default/localassembler.hh:126
Dune::PDELab::LocalAssemblerBase< typename Traits::MatrixBackend, CU, CV > Base
The base class of this local assembler.
Definition: default/localassembler.hh:59
Definition: adaptivity.hh:27
static const bool isNonOverlapping
Definition: default/localassembler.hh:67
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: default/localassembler.hh:169
Definition: assemblerutilities.hh:22
DefaultLocalJacobianAssemblerEngine< DefaultLocalAssembler > LocalJacobianAssemblerEngine
Definition: default/localassembler.hh:83
static bool doLambdaBoundary()
Definition: default/localassembler.hh:187
GO::Traits::TestGridFunctionSpaceConstraints TestGridFunctionSpaceConstraints
The type of the test grid function space constraints.
Definition: assemblerutilities.hh:36
Traits::TestGridFunctionSpace GFSV
Definition: default/localassembler.hh:53
DefaultLocalJacobianApplyAssemblerEngine< DefaultLocalAssembler > LocalJacobianApplyAssemblerEngine
Definition: default/localassembler.hh:84
void setSolution(const Solution &solution_)
Definition: jacobianapplyengine.hh:110
static bool doPatternSkeleton()
Definition: default/localassembler.hh:192
static bool doAlphaVolumePostSkeleton()
Definition: default/localassembler.hh:188
const P & p
Definition: constraints.hh:146
static bool doAlphaVolume()
Query methods for the assembler engines. Theses methods do not belong to the assembler interface...
Definition: default/localassembler.hh:182
GO::Traits::Domain Solution
The type of the domain (solution).
Definition: assemblerutilities.hh:47
Dune::PDELab::LocalFunctionSpace< GFSV, Dune::PDELab::TestSpaceTag > LFSV
Definition: default/localassembler.hh:73
DefaultLocalAssembler(LOP &lop_, const CU &cu_, const CV &cv_, shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
Constructor for non trivial constraints.
Definition: default/localassembler.hh:100
Definition: lfsindexcache.hh:948
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: default/localassembler.hh:140
DefaultLocalPatternAssemblerEngine< DefaultLocalAssembler > LocalPatternAssemblerEngine
Definition: default/localassembler.hh:81
GO::Traits::TestGridFunctionSpace TestGridFunctionSpace
The test grid function space.
Definition: assemblerutilities.hh:29
void postStep()
Definition: default/localassembler.hh:124
RangeField Real
Definition: default/localassembler.hh:50
void setResidual(Residual &residual_)
Definition: default/residualengine.hh:110
void postProcessing(bool v)
Definition: default/localassembler.hh:210
void setSolution(const Solution &solution_)
Definition: default/residualengine.hh:117
static bool doLambdaSkeleton()
Definition: default/localassembler.hh:185
void setResidual(Residual &residual_)
Definition: jacobianapplyengine.hh:103
Traits::Residual::ElementType RangeField
The local operators type for real numbers e.g. time.
Definition: default/localassembler.hh:49
void preStep(Real time_, Real dt_, std::size_t stages_)
Definition: default/localassembler.hh:123