DUNE PDELab (2.8)

residualengine.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
3
5#include <dune/pdelab/gridoperator/common/assemblerutilities.hh>
6#include <dune/pdelab/gridoperator/common/localassemblerenginebase.hh>
7#include <dune/pdelab/constraints/common/constraints.hh>
8#include <dune/pdelab/localoperator/callswitch.hh>
9
10namespace Dune{
11 namespace PDELab{
12
20 template<typename LA>
23 {
24 public:
25
26 template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
27 bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv) const
28 {
29 return false;
30 }
31
33 typedef LA LocalAssembler;
34
36 typedef typename LA::LocalOperator LOP;
37
39 typedef typename LA::Traits::Residual Residual;
40 typedef typename Residual::ElementType ResidualElement;
41
43 typedef typename LA::Traits::Solution Solution;
44 typedef typename Solution::ElementType SolutionElement;
45
47 typedef typename LA::LFSU LFSU;
48 typedef typename LA::LFSUCache LFSUCache;
49 typedef typename LFSU::Traits::GridFunctionSpace GFSU;
50 typedef typename LA::LFSV LFSV;
51 typedef typename LA::LFSVCache LFSVCache;
52 typedef typename LFSV::Traits::GridFunctionSpace GFSV;
53
54 typedef typename Solution::template ConstLocalView<LFSUCache> SolutionView;
55 typedef typename Residual::template LocalView<LFSVCache> ResidualView;
56
64 : local_assembler(local_assembler_),
65 lop(local_assembler_.localOperator()),
66 rl_view(rl,1.0),
67 rn_view(rn,1.0)
68 {}
69
72 bool requireSkeleton() const
73 { return ( local_assembler.doAlphaSkeleton() || local_assembler.doLambdaSkeleton() ); }
74 bool requireSkeletonTwoSided() const
75 { return local_assembler.doSkeletonTwoSided(); }
76 bool requireUVVolume() const
77 { return local_assembler.doAlphaVolume(); }
78 bool requireVVolume() const
79 { return local_assembler.doLambdaVolume(); }
80 bool requireUVSkeleton() const
81 { return local_assembler.doAlphaSkeleton(); }
82 bool requireVSkeleton() const
83 { return local_assembler.doLambdaSkeleton(); }
84 bool requireUVBoundary() const
85 { return local_assembler.doAlphaBoundary(); }
86 bool requireVBoundary() const
87 { return local_assembler.doLambdaBoundary(); }
88 bool requireUVVolumePostSkeleton() const
89 { return local_assembler.doAlphaVolumePostSkeleton(); }
90 bool requireVVolumePostSkeleton() const
91 { return local_assembler.doLambdaVolumePostSkeleton(); }
93
96 {
97 return local_assembler;
98 }
99
101 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
102 {
103 return localAssembler().trialConstraints();
104 }
105
107 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
108 {
109 return localAssembler().testConstraints();
110 }
111
114 void setResidual(Residual & residual_)
115 {
116 global_rl_view.attach(residual_);
117 global_rn_view.attach(residual_);
118 }
119
122 void setSolution(const Solution & solution_)
123 {
124 global_sl_view.attach(solution_);
125 global_sn_view.attach(solution_);
126 }
127
131 template<typename EG, typename LFSUC, typename LFSVC>
132 void onBindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
133 {
134 global_sl_view.bind(lfsu_cache);
135 xl.resize(lfsu_cache.size());
136 }
137
138 template<typename EG, typename LFSVC>
139 void onBindLFSV(const EG & eg, const LFSVC & lfsv_cache)
140 {
141 global_rl_view.bind(lfsv_cache);
142 rl.assign(lfsv_cache.size(),0.0);
143 }
144
145 template<typename IG, typename LFSUC, typename LFSVC>
146 void onBindLFSUVInside(const IG & ig, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
147 {
148 global_sl_view.bind(lfsu_cache);
149 xl.resize(lfsu_cache.size());
150 }
151
152 template<typename IG, typename LFSUC, typename LFSVC>
153 void onBindLFSUVOutside(const IG & ig,
154 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
155 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
156 {
157 global_sn_view.bind(lfsu_n_cache);
158 xn.resize(lfsu_n_cache.size());
159 }
160
161 template<typename IG, typename LFSVC>
162 void onBindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
163 {
164 global_rl_view.bind(lfsv_cache);
165 rl.assign(lfsv_cache.size(),0.0);
166 }
167
168 template<typename IG, typename LFSVC>
169 void onBindLFSVOutside(const IG & ig,
170 const LFSVC & lfsv_s_cache,
171 const LFSVC & lfsv_n_cache)
172 {
173 global_rn_view.bind(lfsv_n_cache);
174 rn.assign(lfsv_n_cache.size(),0.0);
175 }
176
178
182 template<typename EG, typename LFSVC>
183 void onUnbindLFSV(const EG & eg, const LFSVC & lfsv_cache)
184 {
185 global_rl_view.add(rl);
186 global_rl_view.commit();
187 }
188
189 template<typename IG, typename LFSVC>
190 void onUnbindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
191 {
192 global_rl_view.add(rl);
193 global_rl_view.commit();
194 }
195
196 template<typename IG, typename LFSVC>
197 void onUnbindLFSVOutside(const IG & ig,
198 const LFSVC & lfsv_s_cache,
199 const LFSVC & lfsv_n_cache)
200 {
201 global_rn_view.add(rn);
202 global_rn_view.commit();
203 }
205
208 template<typename LFSUC>
209 void loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)
210 {
211 global_sl_view.read(xl);
212 }
213 template<typename LFSUC>
214 void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
215 {
216 global_sn_view.read(xn);
217 }
218 template<typename LFSUC>
219 void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
220 {
221 DUNE_THROW(Dune::NotImplemented,"No coupling lfsu available for ");
222 }
224
227
228 void postAssembly(const GFSU& gfsu, const GFSV& gfsv)
229 {
230 if(local_assembler.doPostProcessing())
231 Dune::PDELab::constrain_residual(local_assembler.testConstraints(),global_rl_view.container());
232
233 }
234
236
239
246 template<typename EG>
247 bool skipEntity(const EG & eg)
248 {
249 return localAssembler().skipEntity(eg);
250 }
251
258 template<typename IG>
259 bool skipIntersection(const IG & ig)
260 {
261 return localAssembler().skipIntersection(ig);
262 }
263
264 template<typename EG, typename LFSUC, typename LFSVC>
265 void assembleUVVolume(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
266 {
267 rl_view.setWeight(local_assembler.weight());
268 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolume>::
269 alpha_volume(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),rl_view);
270 }
271
272 template<typename EG, typename LFSVC>
273 void assembleVVolume(const EG & eg, const LFSVC & lfsv_cache)
274 {
275 rl_view.setWeight(local_assembler.weight());
276 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaVolume>::
277 lambda_volume(lop,eg,lfsv_cache.localFunctionSpace(),rl_view);
278 }
279
280 template<typename IG, typename LFSUC, typename LFSVC>
281 void assembleUVSkeleton(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
282 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
283 {
284 rl_view.setWeight(local_assembler.weight());
285 rn_view.setWeight(local_assembler.weight());
286 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaSkeleton>::
287 alpha_skeleton(lop,ig,
288 lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),
289 lfsu_n_cache.localFunctionSpace(),xn,lfsv_n_cache.localFunctionSpace(),
290 rl_view,rn_view);
291 }
292
293 template<typename IG, typename LFSVC>
294 void assembleVSkeleton(const IG & ig, const LFSVC & lfsv_s_cache, const LFSVC & lfsv_n_cache)
295 {
296 rl_view.setWeight(local_assembler.weight());
297 rn_view.setWeight(local_assembler.weight());
298 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaSkeleton>::
299 lambda_skeleton(lop, ig, lfsv_s_cache.localFunctionSpace(), lfsv_n_cache.localFunctionSpace(), rl_view, rn_view);
300 }
301
302 template<typename IG, typename LFSUC, typename LFSVC>
303 void assembleUVBoundary(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache)
304 {
305 rl_view.setWeight(local_assembler.weight());
306 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaBoundary>::
307 alpha_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),rl_view);
308 }
309
310 template<typename IG, typename LFSVC>
311 void assembleVBoundary(const IG & ig, const LFSVC & lfsv_s_cache)
312 {
313 rl_view.setWeight(local_assembler.weight());
314 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaBoundary>::
315 lambda_boundary(lop,ig,lfsv_s_cache.localFunctionSpace(),rl_view);
316 }
317
318 template<typename IG, typename LFSUC, typename LFSVC>
319 static void assembleUVEnrichedCoupling(const IG & ig,
320 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
321 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache,
322 const LFSUC & lfsu_coupling_cache, const LFSVC & lfsv_coupling_cache)
323 {
324 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
325 }
326
327 template<typename IG, typename LFSVC>
328 static void assembleVEnrichedCoupling(const IG & ig,
329 const LFSVC & lfsv_s_cache,
330 const LFSVC & lfsv_n_cache,
331 const LFSVC & lfsv_coupling_cache)
332 {
333 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
334 }
335
336 template<typename EG, typename LFSUC, typename LFSVC>
337 void assembleUVVolumePostSkeleton(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
338 {
339 rl_view.setWeight(local_assembler.weight());
340 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolumePostSkeleton>::
341 alpha_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),rl_view);
342 }
343
344 template<typename EG, typename LFSVC>
345 void assembleVVolumePostSkeleton(const EG & eg, const LFSVC & lfsv_cache)
346 {
347 rl_view.setWeight(local_assembler.weight());
348 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaVolumePostSkeleton>::
349 lambda_volume_post_skeleton(lop,eg,lfsv_cache.localFunctionSpace(),rl_view);
350 }
351
353
354 private:
357 const LocalAssembler & local_assembler;
358
360 const LOP & lop;
361
363 ResidualView global_rl_view;
364 ResidualView global_rn_view;
365
367 SolutionView global_sl_view;
368 SolutionView global_sn_view;
369
372 typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
373 typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
374
377
379 SolutionVector xl;
381 SolutionVector xn;
383 ResidualVector rl;
385 ResidualVector rn;
391
392 }; // End of class DefaultLocalResidualAssemblerEngine
393
394 }
395}
396#endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
Default exception for dummy implementations.
Definition: exceptions.hh:261
The local assembler engine for DUNE grids which assembles the residual vector.
Definition: residualengine.hh:23
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: residualengine.hh:107
LA::LocalOperator LOP
The type of the local operator.
Definition: residualengine.hh:36
LA LocalAssembler
The type of the wrapping local assembler.
Definition: residualengine.hh:33
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: residualengine.hh:132
LA::Traits::Residual Residual
The type of the residual vector.
Definition: residualengine.hh:39
LA::LFSU LFSU
The local function spaces.
Definition: residualengine.hh:47
void setResidual(Residual &residual_)
Definition: residualengine.hh:114
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: residualengine.hh:183
bool requireSkeleton() const
Definition: residualengine.hh:72
bool skipIntersection(const IG &ig)
Definition: residualengine.hh:259
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: residualengine.hh:228
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: residualengine.hh:95
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: residualengine.hh:101
LA::Traits::Solution Solution
The type of the solution vector.
Definition: residualengine.hh:43
bool skipEntity(const EG &eg)
Definition: residualengine.hh:247
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: residualengine.hh:209
void setSolution(const Solution &solution_)
Definition: residualengine.hh:122
DefaultLocalResidualAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: residualengine.hh:63
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:22
void assign(size_type size, const T &value)
Resize the container to size and assign the passed value to all entries.
Definition: localvector.hh:330
WeightedVectorAccumulationView< LocalVector > WeightedAccumulationView
An accumulate-only view of this container that automatically applies a weight to all contributions.
Definition: localvector.hh:211
void resize(size_type size)
Resize the container.
Definition: localvector.hh:324
void setWeight(weight_type weight)
Resets the weighting coefficient of the view.
Definition: localvector.hh:72
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:904
Dune namespace.
Definition: alignedallocator.hh:11
Definition: localfunctionspacetags.hh:54
Definition: localfunctionspacetags.hh:48
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)