DUNE PDELab (git)

localassembler.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_LOCALASSEMBLER_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_LOCALASSEMBLER_HH
3
4#include <dune/typetree/typetree.hh>
5
6#include <dune/pdelab/gridoperator/default/residualengine.hh>
7#include <dune/pdelab/gridoperator/default/patternengine.hh>
8#include <dune/pdelab/gridoperator/default/jacobianengine.hh>
9#include <dune/pdelab/gridoperator/default/jacobianapplyengine.hh>
10#include <dune/pdelab/gridoperator/common/assemblerutilities.hh>
11#include <dune/pdelab/gridfunctionspace/lfsindexcache.hh>
12
13namespace 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 public:
37
40
42 typedef typename Traits::Residual::ElementType RangeField;
43 typedef RangeField Real;
44
45 typedef typename Traits::TrialGridFunctionSpace GFSU;
46 typedef typename Traits::TestGridFunctionSpace GFSV;
47
50
53
55 typedef LOP LocalOperator;
56
57 static const bool isNonOverlapping = nonoverlapping_mode;
58
61 // Types of local function spaces
64 typedef LFSIndexCache<LFSU,CU> LFSUCache;
65 typedef LFSIndexCache<LFSV,CV> LFSVCache;
66
68
75
76 // friend declarations such that engines are able to call scatter_jacobian() and add_entry() from base class
80
82 DefaultLocalAssembler (LOP & lop, std::shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
83 : lop_(lop), weight_(1.0), doPreProcessing_(true), doPostProcessing_(true),
84 pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this)
85 , jacobian_apply_engine(*this)
86 , _reconstruct_border_entries(isNonOverlapping)
87 {}
88
90 DefaultLocalAssembler (LOP & lop, const CU& cu, const CV& cv,
91 std::shared_ptr<typename GO::BorderDOFExchanger> border_dof_exchanger)
92 : Base(cu, cv),
93 lop_(lop), weight_(1.0), doPreProcessing_(true), doPostProcessing_(true),
94 pattern_engine(*this,border_dof_exchanger), residual_engine(*this), jacobian_engine(*this)
95 , jacobian_apply_engine(*this)
96 , _reconstruct_border_entries(isNonOverlapping)
97 {}
98
101 {
102 return lop_;
103 }
104
106 const LOP &localOperator() const
107 {
108 return lop_;
109 }
110
117 template<class EG>
118 bool skipEntity(const EG & eg) const
119 {
120 static_assert(
121 models<Impl::HasDoSkipEntity, LOP>(),
122 "Your local operator does not provide the 'doSkipEntity' flag. "
123 "If you are porting a previous implementation, set the flag to "
124 "false");
125 return Dune::PDELab::LocalAssemblerCallSwitch<
126 LOP, LOP::doSkipEntity>::skip_entity(lop_, eg);
127 }
128
135 template<class IG>
136 bool skipIntersection(const IG & ig) const
137 {
138 static_assert(models<Impl::HasDoSkipIntersection, LOP>(),
139 "Your local operator does not provide the "
140 "'doSkipIntersection' flag. "
141 "If you are porting a previous implementation, set the "
142 "flag to false");
143 return Dune::PDELab::LocalAssemblerCallSwitch<
144 LOP, LOP::doSkipIntersection>::skip_intersection(lop_, ig);
145 }
146
150 void setTime(Real time_)
151 {
152 lop_.setTime(time_);
153 }
154
157 {
158 return weight_;
159 }
160
163 weight_ = weight;
164 }
165
168 void preStage (Real time_, int r_) { lop_.preStage(time_,r_); }
169 void preStep (Real time_, Real dt_, std::size_t stages_){ lop_.preStep(time_,dt_,stages_); }
170 void postStep (){ lop_.postStep(); }
171 void postStage (){ lop_.postStage(); }
172 Real suggestTimestep (Real dt) const{return lop_.suggestTimestep(dt); }
174
175 bool reconstructBorderEntries() const
176 {
177 return _reconstruct_border_entries;
178 }
179
182
186 (typename Traits::MatrixPattern & p)
187 {
188 pattern_engine.setPattern(p);
189 return pattern_engine;
190 }
191
195 (typename Traits::Residual & r, const typename Traits::Solution & x)
196 {
197 residual_engine.setResidual(r);
198 residual_engine.setSolution(x);
199 return residual_engine;
200 }
201
205 (typename Traits::Jacobian & a, const typename Traits::Solution & x)
206 {
207 jacobian_engine.setJacobian(a);
208 jacobian_engine.setSolution(x);
209 return jacobian_engine;
210 }
211
215 (const typename Traits::Domain & update, typename Traits::Range & result)
216 {
217 jacobian_apply_engine.setUpdate(update);
218 jacobian_apply_engine.setResult(result);
219 return jacobian_apply_engine;
220 }
221
225 (const typename Traits::Domain & solution, const typename Traits::Domain & update, typename Traits::Range & result)
226 {
227 jacobian_apply_engine.setSolution(solution);
228 jacobian_apply_engine.setUpdate(update);
229 jacobian_apply_engine.setResult(result);
230 return jacobian_apply_engine;
231 }
232
234
239 static constexpr bool doAlphaVolume() { return LOP::doAlphaVolume; }
240 static constexpr bool doLambdaVolume() { return LOP::doLambdaVolume; }
241 static constexpr bool doAlphaSkeleton() { return LOP::doAlphaSkeleton; }
242 static constexpr bool doLambdaSkeleton() { return LOP::doLambdaSkeleton; }
243 static constexpr bool doAlphaBoundary() { return LOP::doAlphaBoundary; }
244 static constexpr bool doLambdaBoundary() { return LOP::doLambdaBoundary; }
245 static constexpr bool doAlphaVolumePostSkeleton() { return LOP::doAlphaVolumePostSkeleton; }
246 static constexpr bool doLambdaVolumePostSkeleton() { return LOP::doLambdaVolumePostSkeleton; }
247 static constexpr bool doSkeletonTwoSided() { return LOP::doSkeletonTwoSided; }
248 static constexpr bool doPatternVolume() { return LOP::doPatternVolume; }
249 static constexpr bool doPatternSkeleton() { return LOP::doPatternSkeleton; }
250 static constexpr bool doPatternBoundary() { return LOP::doPatternBoundary; }
251 static constexpr bool doPatternVolumePostSkeleton() { return LOP::doPatternVolumePostSkeleton; }
252 static constexpr bool isLinear() { return LOP::isLinear;}
254
256
259 bool doPreProcessing() const { return doPreProcessing_; }
260
265 void preProcessing(bool v)
266 {
267 doPreProcessing_ = v;
268 }
269
271
274 bool doPostProcessing() const { return doPostProcessing_; }
275
280 void postProcessing(bool v)
281 {
282 doPostProcessing_ = v;
283 }
284
285 private:
286
288 LOP & lop_;
289
291 RangeField weight_;
292
295 bool doPreProcessing_;
296
299 bool doPostProcessing_;
300
303 LocalPatternAssemblerEngine pattern_engine;
304 LocalResidualAssemblerEngine residual_engine;
305 LocalJacobianAssemblerEngine jacobian_engine;
306 LocalJacobianApplyAssemblerEngine jacobian_apply_engine;
308
309 bool _reconstruct_border_entries;
310 };
311
312 } // end namespace PDELab
313} // end namespace Dune
314#endif
The local assembler for DUNE grids.
Definition: localassembler.hh:35
Traits::Residual::ElementType RangeField
The local operators type for real numbers e.g. time.
Definition: localassembler.hh:42
RangeField weight() const
Obtain the weight that was set last.
Definition: localassembler.hh:156
DefaultLocalPatternAssemblerEngine< DefaultLocalAssembler > LocalPatternAssemblerEngine
Definition: localassembler.hh:71
Dune::PDELab::LocalFunctionSpace< GFSU, Dune::PDELab::TrialSpaceTag > LFSU
Definition: localassembler.hh:62
static constexpr bool doAlphaVolume()
Query methods for the assembler engines. Theses methods do not belong to the assembler interface,...
Definition: localassembler.hh:239
void setTime(Real time_)
Definition: localassembler.hh:150
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: localassembler.hh:186
bool doPreProcessing() const
Query whether to do preprocessing in the engines.
Definition: localassembler.hh:259
Dune::PDELab::LocalAssemblerTraits< GO > Traits
The traits class.
Definition: localassembler.hh:39
void preStage(Real time_, int r_)
Definition: localassembler.hh:168
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(const typename Traits::Domain &update, typename Traits::Range &result)
Definition: localassembler.hh:215
const LOP & localOperator() const
get a reference to the local operator
Definition: localassembler.hh:106
LOP LocalOperator
The local operator.
Definition: localassembler.hh:55
bool skipIntersection(const IG &ig) const
Definition: localassembler.hh:136
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: localassembler.hh:195
bool skipEntity(const EG &eg) const
Definition: localassembler.hh:118
void preProcessing(bool v)
Definition: localassembler.hh:265
LOP & localOperator()
get a reference to the local operator
Definition: localassembler.hh:100
void postProcessing(bool v)
Definition: localassembler.hh:280
void setWeight(RangeField weight)
Notifies the assembler about the current weight of assembling.
Definition: localassembler.hh:162
DefaultLocalAssembler(LOP &lop, std::shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
Constructor with empty constraints.
Definition: localassembler.hh:82
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: localassembler.hh:205
Dune::PDELab::LocalAssemblerBase< typename Traits::MatrixBackend, CU, CV > Base
The base class of this local assembler.
Definition: localassembler.hh:52
DefaultLocalAssembler(LOP &lop, const CU &cu, const CV &cv, std::shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
Constructor for non trivial constraints.
Definition: localassembler.hh:90
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(const typename Traits::Domain &solution, const typename Traits::Domain &update, typename Traits::Range &result)
Definition: localassembler.hh:225
bool doPostProcessing() const
Query whether to do postprocessing in the engines.
Definition: localassembler.hh:274
void setUpdate(const Domain &update_)
Definition: jacobianapplyengine.hh:117
void setResult(Range &result_)
Definition: jacobianapplyengine.hh:125
void setSolution(const Domain &solution_)
Definition: jacobianapplyengine.hh:109
void setSolution(const Solution &solution_)
Definition: jacobianengine.hh:120
void setJacobian(Jacobian &jacobian_)
Definition: jacobianengine.hh:110
void setPattern(Pattern &pattern_)
Definition: patternengine.hh:92
void setResidual(Residual &residual_)
Definition: residualengine.hh:114
void setSolution(const Solution &solution_)
Definition: residualengine.hh:122
Base class for local assembler.
Definition: assemblerutilities.hh:217
Create a local function space from a global function space.
Definition: localfunctionspace.hh:754
Dune namespace.
Definition: alignedallocator.hh:13
Definition: assemblerutilities.hh:51
GO::Traits::Range Residual
The type of the range (residual).
Definition: assemblerutilities.hh:88
MatrixBackend::template Pattern< Jacobian, TestGridFunctionSpace, TrialGridFunctionSpace > MatrixPattern
The matrix pattern.
Definition: assemblerutilities.hh:102
GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace
The trial grid function space.
Definition: assemblerutilities.hh:54
GO::Traits::Jacobian Jacobian
The type of the jacobian.
Definition: assemblerutilities.hh:95
GO::Traits::TestGridFunctionSpace TestGridFunctionSpace
The test grid function space.
Definition: assemblerutilities.hh:57
GO::Traits::Domain Solution
The type of the domain (solution).
Definition: assemblerutilities.hh:78
GO::Traits::Range Range
The type of the range (residual).
Definition: assemblerutilities.hh:85
GO::Traits::Domain Domain
The type of the domain (solution).
Definition: assemblerutilities.hh:75
GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints
The type of the trial grid function space constraints.
Definition: assemblerutilities.hh:61
GO::Traits::TestGridFunctionSpaceConstraints TestGridFunctionSpaceConstraints
The type of the test grid function space constraints.
Definition: assemblerutilities.hh:64
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)