dune-pdelab  2.5-dev
default/localassembler.hh
Go to the documentation of this file.
1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_LOCALASSEMBLER_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_LOCALASSEMBLER_HH
3 
4 #include <dune/typetree/typetree.hh>
5 
13 
14 namespace Dune{
15  namespace PDELab{
16 
31  template<typename GO, typename LOP, bool nonoverlapping_mode = false>
33  public Dune::PDELab::LocalAssemblerBase<typename GO::Traits::MatrixBackend,
34  typename GO::Traits::TrialGridFunctionSpaceConstraints,
35  typename GO::Traits::TestGridFunctionSpaceConstraints>
36  {
37  public:
38 
41 
43  typedef typename Traits::Residual::ElementType RangeField;
44  typedef RangeField Real;
45 
48 
51 
54 
56  typedef LOP LocalOperator;
57 
58  static const bool isNonOverlapping = nonoverlapping_mode;
59 
62  // Types of local function spaces
67 
69 
77 
78  // friend declarations such that engines are able to call scatter_jacobian() and add_entry() from base class
82 
84  DefaultLocalAssembler (LOP & lop, shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
85  : lop_(lop), weight_(1.0), doPreProcessing_(true), doPostProcessing_(true),
86  pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this)
87  , jacobian_apply_engine(*this)
88  , nonlinear_jacobian_apply_engine(*this)
89  , _reconstruct_border_entries(isNonOverlapping)
90  {}
91 
93  DefaultLocalAssembler (LOP & lop, const CU& cu, const CV& cv,
94  shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
95  : Base(cu, cv),
96  lop_(lop), weight_(1.0), doPreProcessing_(true), doPostProcessing_(true),
97  pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this)
98  , jacobian_apply_engine(*this)
99  , nonlinear_jacobian_apply_engine(*this)
100  , _reconstruct_border_entries(isNonOverlapping)
101  {}
102 
105  {
106  return lop_;
107  }
108 
110  const LOP &localOperator() const
111  {
112  return lop_;
113  }
114 
118  void setTime(Real time_)
119  {
120  lop_.setTime(time_);
121  }
122 
125  {
126  return weight_;
127  }
128 
131  weight_ = weight;
132  }
133 
136  void preStage (Real time_, int r_) { lop_.preStage(time_,r_); }
137  void preStep (Real time_, Real dt_, std::size_t stages_){ lop_.preStep(time_,dt_,stages_); }
138  void postStep (){ lop_.postStep(); }
139  void postStage (){ lop_.postStage(); }
140  Real suggestTimestep (Real dt) const{return lop_.suggestTimestep(dt); }
142 
144  {
145  return _reconstruct_border_entries;
146  }
147 
150 
155  {
156  pattern_engine.setPattern(p);
157  return pattern_engine;
158  }
159 
163  (typename Traits::Residual & r, const typename Traits::Solution & x)
164  {
165  residual_engine.setResidual(r);
166  residual_engine.setSolution(x);
167  return residual_engine;
168  }
169 
173  (typename Traits::Jacobian & a, const typename Traits::Solution & x)
174  {
175  jacobian_engine.setJacobian(a);
176  jacobian_engine.setSolution(x);
177  return jacobian_engine;
178  }
179 
183  (typename Traits::Residual & r, const typename Traits::Solution & x)
184  {
185  jacobian_apply_engine.setResidual(r);
186  jacobian_apply_engine.setSolution(x);
187  return jacobian_apply_engine;
188  }
189 
193  (typename Traits::Residual & r, const typename Traits::Solution & x, const typename Traits::Solution & z)
194  {
195  nonlinear_jacobian_apply_engine.setResidual(r);
196  nonlinear_jacobian_apply_engine.setSolution(x);
197  nonlinear_jacobian_apply_engine.setUpdate(z);
198  return nonlinear_jacobian_apply_engine;
199  }
200 
202 
207  static bool doAlphaVolume() { return LOP::doAlphaVolume; }
208  static bool doLambdaVolume() { return LOP::doLambdaVolume; }
209  static bool doAlphaSkeleton() { return LOP::doAlphaSkeleton; }
210  static bool doLambdaSkeleton() { return LOP::doLambdaSkeleton; }
211  static bool doAlphaBoundary() { return LOP::doAlphaBoundary; }
212  static bool doLambdaBoundary() { return LOP::doLambdaBoundary; }
213  static bool doAlphaVolumePostSkeleton() { return LOP::doAlphaVolumePostSkeleton; }
214  static bool doLambdaVolumePostSkeleton() { return LOP::doLambdaVolumePostSkeleton; }
215  static bool doSkeletonTwoSided() { return LOP::doSkeletonTwoSided; }
216  static bool doPatternVolume() { return LOP::doPatternVolume; }
217  static bool doPatternSkeleton() { return LOP::doPatternSkeleton; }
218  static bool doPatternBoundary() { return LOP::doPatternBoundary; }
219  static bool doPatternVolumePostSkeleton() { return LOP::doPatternVolumePostSkeleton; }
221 
223 
226  bool doPreProcessing() const { return doPreProcessing_; }
227 
232  void preProcessing(bool v)
233  {
234  doPreProcessing_ = v;
235  }
236 
238 
241  bool doPostProcessing() const { return doPostProcessing_; }
242 
247  void postProcessing(bool v)
248  {
249  doPostProcessing_ = v;
250  }
251 
252  private:
253 
255  LOP & lop_;
256 
258  RangeField weight_;
259 
262  bool doPreProcessing_;
263 
266  bool doPostProcessing_;
267 
270  LocalPatternAssemblerEngine pattern_engine;
271  LocalResidualAssemblerEngine residual_engine;
272  LocalJacobianAssemblerEngine jacobian_engine;
273  LocalJacobianApplyAssemblerEngine jacobian_apply_engine;
274  LocalNonlinearJacobianApplyAssemblerEngine nonlinear_jacobian_apply_engine;
276 
277  bool _reconstruct_border_entries;
278 
279  };
280 
281  }
282 }
283 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_LOCALASSEMBLER_HH
Dune::PDELab::DefaultLocalAssembler::doPatternVolumePostSkeleton
static bool doPatternVolumePostSkeleton()
Definition: default/localassembler.hh:219
assemblerutilities.hh
Dune::PDELab::DefaultLocalAssembler::doLambdaBoundary
static bool doLambdaBoundary()
Definition: default/localassembler.hh:212
Dune::PDELab::DefaultLocalAssembler::Base
Dune::PDELab::LocalAssemblerBase< typename Traits::MatrixBackend, CU, CV > Base
The base class of this local assembler.
Definition: default/localassembler.hh:53
Dune::PDELab::DefaultLocalAssembler::doAlphaSkeleton
static bool doAlphaSkeleton()
Definition: default/localassembler.hh:209
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine< DefaultLocalAssembler >
Dune::PDELab::DefaultLocalAssembler::CU
Traits::TrialGridFunctionSpaceConstraints CU
Definition: default/localassembler.hh:49
Dune::PDELab::LocalAssemblerBase
Base class for local assembler.
Definition: assemblerutilities.hh:182
Dune::PDELab::DefaultLocalAssembler::reconstructBorderEntries
bool reconstructBorderEntries() const
Definition: default/localassembler.hh:143
Dune::PDELab::DefaultLocalAssembler::preStage
void preStage(Real time_, int r_)
Definition: default/localassembler.hh:136
lfsindexcache.hh
Dune::PDELab::LocalAssemblerTraits::Solution
GO::Traits::Domain Solution
The type of the domain (solution).
Definition: assemblerutilities.hh:47
Dune::PDELab::DefaultLocalAssembler::postStage
void postStage()
Definition: default/localassembler.hh:139
Dune::PDELab::DefaultLocalAssembler::localPatternAssemblerEngine
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: default/localassembler.hh:154
Dune::PDELab::DefaultLocalAssembler::weight
RangeField weight() const
Obtain the weight that was set last.
Definition: default/localassembler.hh:124
Dune::PDELab::DefaultLocalResidualAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: default/residualengine.hh:122
Dune::PDELab::DefaultLocalAssembler::LFSUCache
LFSIndexCache< LFSU, CU > LFSUCache
Definition: default/localassembler.hh:65
Dune::PDELab::DefaultLocalAssembler::preStep
void preStep(Real time_, Real dt_, std::size_t stages_)
Definition: default/localassembler.hh:137
Dune::PDELab::DefaultLocalAssembler::LocalOperator
LOP LocalOperator
The local operator.
Definition: default/localassembler.hh:56
Dune
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Dune::PDELab::DefaultLocalAssembler::localOperator
const LOP & localOperator() const
get a reference to the local operator
Definition: default/localassembler.hh:110
Dune::PDELab::DefaultLocalAssembler::LocalJacobianAssemblerEngine
DefaultLocalJacobianAssemblerEngine< DefaultLocalAssembler > LocalJacobianAssemblerEngine
Definition: default/localassembler.hh:74
Dune::PDELab::DefaultLocalAssembler::setWeight
void setWeight(RangeField weight)
Notifies the assembler about the current weight of assembling.
Definition: default/localassembler.hh:130
Dune::PDELab::DefaultLocalAssembler::LFSVCache
LFSIndexCache< LFSV, CV > LFSVCache
Definition: default/localassembler.hh:66
Dune::PDELab::DefaultLocalAssembler::RangeField
Traits::Residual::ElementType RangeField
The local operators type for real numbers e.g. time.
Definition: default/localassembler.hh:43
Dune::PDELab::DefaultLocalAssembler::doAlphaBoundary
static bool doAlphaBoundary()
Definition: default/localassembler.hh:211
Dune::PDELab::DefaultLocalJacobianApplyAssemblerEngine::setResidual
void setResidual(Residual &residual_)
Definition: default/jacobianapplyengine.hh:106
Dune::PDELab::DefaultLocalAssembler::doLambdaVolumePostSkeleton
static bool doLambdaVolumePostSkeleton()
Definition: default/localassembler.hh:214
Dune::PDELab::DefaultLocalJacobianApplyAssemblerEngine< DefaultLocalAssembler >
residualengine.hh
Dune::PDELab::DefaultLocalAssembler::doLambdaSkeleton
static bool doLambdaSkeleton()
Definition: default/localassembler.hh:210
Dune::PDELab::DefaultLocalAssembler::suggestTimestep
Real suggestTimestep(Real dt) const
Definition: default/localassembler.hh:140
Dune::PDELab::LocalAssemblerTraits
Definition: assemblerutilities.hh:22
Dune::PDELab::DefaultLocalAssembler::DefaultLocalAssembler
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:93
Dune::PDELab::DefaultLocalAssembler::DefaultLocalAssembler
DefaultLocalAssembler(LOP &lop, shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
Constructor with empty constraints.
Definition: default/localassembler.hh:84
Dune::PDELab::DefaultLocalAssembler::postStep
void postStep()
Definition: default/localassembler.hh:138
Dune::PDELab::DefaultLocalAssembler::doPatternVolume
static bool doPatternVolume()
Definition: default/localassembler.hh:216
Dune::PDELab::LocalAssemblerTraits::TrialGridFunctionSpaceConstraints
GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints
The type of the trial grid function space constraints.
Definition: assemblerutilities.hh:33
Dune::PDELab::DefaultLocalAssembler::doSkeletonTwoSided
static bool doSkeletonTwoSided()
Definition: default/localassembler.hh:215
Dune::PDELab::DefaultLocalAssembler::CV
Traits::TestGridFunctionSpaceConstraints CV
Definition: default/localassembler.hh:50
jacobianengine.hh
Dune::PDELab::DefaultLocalAssembler::doLambdaVolume
static bool doLambdaVolume()
Definition: default/localassembler.hh:208
Dune::PDELab::DefaultLocalAssembler::localResidualAssemblerEngine
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: default/localassembler.hh:163
Dune::PDELab::DefaultLocalAssembler::GFSU
Traits::TrialGridFunctionSpace GFSU
Definition: default/localassembler.hh:46
Dune::PDELab::DefaultLocalAssembler::Traits
Dune::PDELab::LocalAssemblerTraits< GO > Traits
The traits class.
Definition: default/localassembler.hh:40
Dune::PDELab::DefaultLocalAssembler::setTime
void setTime(Real time_)
Definition: default/localassembler.hh:118
Dune::PDELab::DefaultLocalAssembler::doAlphaVolumePostSkeleton
static bool doAlphaVolumePostSkeleton()
Definition: default/localassembler.hh:213
Dune::PDELab::DefaultLocalAssembler::doPatternBoundary
static bool doPatternBoundary()
Definition: default/localassembler.hh:218
Dune::PDELab::LocalAssemblerTraits::TrialGridFunctionSpace
GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace
The trial grid function space.
Definition: assemblerutilities.hh:26
Dune::PDELab::DefaultLocalPatternAssemblerEngine< DefaultLocalAssembler >
patternengine.hh
Dune::PDELab::DefaultLocalAssembler::GFSV
Traits::TestGridFunctionSpace GFSV
Definition: default/localassembler.hh:47
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: default/nonlinearjacobianapplyengine.hh:113
Dune::PDELab::DefaultLocalAssembler::LFSU
Dune::PDELab::LocalFunctionSpace< GFSU, Dune::PDELab::TrialSpaceTag > LFSU
Definition: default/localassembler.hh:63
Dune::PDELab::LocalAssemblerTraits::MatrixPattern
MatrixBackend::template Pattern< Jacobian, TestGridFunctionSpace, TrialGridFunctionSpace > MatrixPattern
The matrix pattern.
Definition: assemblerutilities.hh:68
Dune::PDELab::DefaultLocalAssembler
The local assembler for DUNE grids.
Definition: default/localassembler.hh:32
Dune::PDELab::LocalAssemblerTraits::Residual
GO::Traits::Range Residual
The type of the range (residual).
Definition: assemblerutilities.hh:54
Dune::PDELab::DefaultLocalAssembler::localJacobianAssemblerEngine
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: default/localassembler.hh:173
Dune::PDELab::DefaultLocalAssembler::Real
RangeField Real
Definition: default/localassembler.hh:44
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setUpdate
void setUpdate(const Solution &update_)
Definition: default/nonlinearjacobianapplyengine.hh:121
Dune::PDELab::DefaultLocalPatternAssemblerEngine::setPattern
void setPattern(Pattern &pattern_)
Definition: default/patternengine.hh:92
Dune::PDELab::DefaultLocalAssembler::LocalPatternAssemblerEngine
DefaultLocalPatternAssemblerEngine< DefaultLocalAssembler > LocalPatternAssemblerEngine
Definition: default/localassembler.hh:72
Dune::PDELab::LocalAssemblerTraits::TestGridFunctionSpaceConstraints
GO::Traits::TestGridFunctionSpaceConstraints TestGridFunctionSpaceConstraints
The type of the test grid function space constraints.
Definition: assemblerutilities.hh:36
Dune::PDELab::DefaultLocalJacobianApplyAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: default/jacobianapplyengine.hh:114
Dune::PDELab::DefaultLocalAssembler::localOperator
LOP & localOperator()
get a reference to the local operator
Definition: default/localassembler.hh:104
Dune::PDELab::DefaultLocalAssembler::doAlphaVolume
static bool doAlphaVolume()
Query methods for the assembler engines. Theses methods do not belong to the assembler interface,...
Definition: default/localassembler.hh:207
Dune::PDELab::DefaultLocalAssembler::doPostProcessing
bool doPostProcessing() const
Query whether to do postprocessing in the engines.
Definition: default/localassembler.hh:241
Dune::PDELab::DefaultLocalAssembler::localNonlinearJacobianApplyAssemblerEngine
LocalNonlinearJacobianApplyAssemblerEngine & localNonlinearJacobianApplyAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x, const typename Traits::Solution &z)
Definition: default/localassembler.hh:193
Dune::PDELab::DefaultLocalAssembler::LocalJacobianApplyAssemblerEngine
DefaultLocalJacobianApplyAssemblerEngine< DefaultLocalAssembler > LocalJacobianApplyAssemblerEngine
Definition: default/localassembler.hh:75
Dune::PDELab::LocalAssemblerTraits::TestGridFunctionSpace
GO::Traits::TestGridFunctionSpace TestGridFunctionSpace
The test grid function space.
Definition: assemblerutilities.hh:29
Dune::PDELab::DefaultLocalAssembler::LocalResidualAssemblerEngine
DefaultLocalResidualAssemblerEngine< DefaultLocalAssembler > LocalResidualAssemblerEngine
Definition: default/localassembler.hh:73
Dune::PDELab::DefaultLocalJacobianAssemblerEngine< DefaultLocalAssembler >
jacobianapplyengine.hh
Dune::PDELab::DefaultLocalResidualAssemblerEngine::setResidual
void setResidual(Residual &residual_)
Definition: default/residualengine.hh:114
Dune::PDELab::DefaultLocalAssembler::LFSV
Dune::PDELab::LocalFunctionSpace< GFSV, Dune::PDELab::TestSpaceTag > LFSV
Definition: default/localassembler.hh:64
Dune::PDELab::DefaultLocalAssembler::doPatternSkeleton
static bool doPatternSkeleton()
Definition: default/localassembler.hh:217
Dune::PDELab::DefaultLocalAssembler::postProcessing
void postProcessing(bool v)
Definition: default/localassembler.hh:247
Dune::PDELab::LocalFunctionSpace
Create a local function space from a global function space.
Definition: localfunctionspace.hh:698
Dune::PDELab::DefaultLocalJacobianAssemblerEngine::setSolution
void setSolution(const Solution &solution_)
Definition: default/jacobianengine.hh:120
Dune::PDELab::DefaultLocalAssembler::preProcessing
void preProcessing(bool v)
Definition: default/localassembler.hh:232
Dune::PDELab::DefaultLocalJacobianAssemblerEngine::setJacobian
void setJacobian(Jacobian &jacobian_)
Definition: default/jacobianengine.hh:110
Dune::PDELab::DefaultLocalNonlinearJacobianApplyAssemblerEngine::setResidual
void setResidual(Residual &residual_)
Definition: default/nonlinearjacobianapplyengine.hh:105
Dune::PDELab::DefaultLocalAssembler::doPreProcessing
bool doPreProcessing() const
Query whether to do preprocessing in the engines.
Definition: default/localassembler.hh:226
nonlinearjacobianapplyengine.hh
Dune::PDELab::DefaultLocalAssembler::localJacobianApplyAssemblerEngine
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: default/localassembler.hh:183
Dune::PDELab::DefaultLocalResidualAssemblerEngine< DefaultLocalAssembler >
Dune::PDELab::LocalAssemblerTraits::Jacobian
GO::Traits::Jacobian Jacobian
The type of the jacobian.
Definition: assemblerutilities.hh:61
Dune::PDELab::DefaultLocalAssembler::LocalNonlinearJacobianApplyAssemblerEngine
DefaultLocalNonlinearJacobianApplyAssemblerEngine< DefaultLocalAssembler > LocalNonlinearJacobianApplyAssemblerEngine
Definition: default/localassembler.hh:76
Dune::PDELab::DefaultLocalAssembler::isNonOverlapping
static const bool isNonOverlapping
Definition: default/localassembler.hh:58
p
const P & p
Definition: constraints.hh:147
Dune::PDELab::LFSIndexCache
Definition: lfsindexcache.hh:977