DUNE PDELab (git)

residualengine.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_FASTDG_RESIDUALENGINE_HH
2#define DUNE_PDELAB_GRIDOPERATOR_FASTDG_RESIDUALENGINE_HH
3
5#include <dune/pdelab/gridoperator/common/gridoperatorutilities.hh>
6#include <dune/pdelab/gridoperator/common/assemblerutilities.hh>
7#include <dune/pdelab/gridoperator/common/localassemblerenginebase.hh>
8#include <dune/pdelab/constraints/common/constraints.hh>
9#include <dune/pdelab/localoperator/callswitch.hh>
10
11namespace Dune{
12 namespace PDELab{
13
21 template<typename LA>
24 {
25 public:
26
27 template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
28 bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv) const
29 {
30 return false;
31 }
32
34 typedef LA LocalAssembler;
35
37 typedef typename LA::LocalOperator LOP;
38
40 typedef typename LA::Traits::Residual Residual;
41 typedef typename Residual::ElementType ResidualElement;
42
44 typedef typename LA::Traits::Solution Solution;
45 typedef typename Solution::ElementType SolutionElement;
46
48 typedef typename LA::LFSU LFSU;
49 typedef typename LA::LFSUCache LFSUCache;
50 typedef typename LFSU::Traits::GridFunctionSpace GFSU;
51 typedef typename LA::LFSV LFSV;
52 typedef typename LA::LFSVCache LFSVCache;
53 typedef typename LFSV::Traits::GridFunctionSpace GFSV;
54
55 typedef typename Solution::template ConstAliasedLocalView<LFSUCache> SolutionView;
56 typedef typename Residual::template AliasedLocalView<LFSVCache> ResidualView;
57
65 : local_assembler(local_assembler_),
66 lop(local_assembler_.localOperator())
67 //rl_view(rl,1.0),
68 //rn_view(rn,1.0)
69 {}
70
72
87 local_assembler(other.local_assembler), lop(other.lop),
88 global_rl_view(other.global_rl_view),
89 global_rn_view(other.global_rn_view),
90 global_sl_view(other.global_sl_view),
91 global_sn_view(other.global_sn_view)
92 //rl_view(rl,1.0),
93 //rn_view(rn,1.0)
94 { }
95
98 bool requireSkeleton() const
99 { return ( local_assembler.doAlphaSkeleton() || local_assembler.doLambdaSkeleton() ); }
100 bool requireSkeletonTwoSided() const
101 { return local_assembler.doSkeletonTwoSided(); }
102 bool requireUVVolume() const
103 { return local_assembler.doAlphaVolume(); }
104 bool requireVVolume() const
105 { return local_assembler.doLambdaVolume(); }
106 bool requireUVSkeleton() const
107 { return local_assembler.doAlphaSkeleton(); }
108 bool requireVSkeleton() const
109 { return local_assembler.doLambdaSkeleton(); }
110 bool requireUVBoundary() const
111 { return local_assembler.doAlphaBoundary(); }
112 bool requireVBoundary() const
113 { return local_assembler.doLambdaBoundary(); }
114 bool requireUVVolumePostSkeleton() const
115 { return local_assembler.doAlphaVolumePostSkeleton(); }
116 bool requireVVolumePostSkeleton() const
117 { return local_assembler.doLambdaVolumePostSkeleton(); }
119
122 {
123 return local_assembler;
124 }
125
127 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
128 {
129 return localAssembler().trialConstraints();
130 }
131
133 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
134 {
135 return localAssembler().testConstraints();
136 }
137
140 void setResidual(Residual & residual_)
141 {
142 global_rl_view.attach(residual_);
143 global_rn_view.attach(residual_);
144 }
145
148 void setSolution(const Solution & solution_)
149 {
150 global_sl_view.attach(solution_);
151 global_sn_view.attach(solution_);
152 }
153
157 template<typename EG, typename LFSUC, typename LFSVC>
158 void onBindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
159 {
160 global_sl_view.bind(lfsu_cache);
161 //xl.resize(lfsu_cache.size());
162 }
163
164 template<typename EG, typename LFSVC>
165 void onBindLFSV(const EG & eg, const LFSVC & lfsv_cache)
166 {
167 global_rl_view.bind(lfsv_cache);
168 //rl.assign(lfsv_cache.size(),0.0);
169 }
170
171 template<typename IG, typename LFSUC, typename LFSVC>
172 void onBindLFSUVInside(const IG & ig, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
173 {
174 global_sl_view.bind(lfsu_cache);
175 //xl.resize(lfsu_cache.size());
176 }
177
178 template<typename IG, typename LFSUC, typename LFSVC>
179 void onBindLFSUVOutside(const IG & ig,
180 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
181 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
182 {
183 global_sn_view.bind(lfsu_n_cache);
184 //xn.resize(lfsu_n_cache.size());
185 }
186
187 template<typename IG, typename LFSVC>
188 void onBindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
189 {
190 global_rl_view.bind(lfsv_cache);
191 //rl.assign(lfsv_cache.size(),0.0);
192 }
193
194 template<typename IG, typename LFSVC>
195 void onBindLFSVOutside(const IG & ig,
196 const LFSVC & lfsv_s_cache,
197 const LFSVC & lfsv_n_cache)
198 {
199 global_rn_view.bind(lfsv_n_cache);
200 //rn.assign(lfsv_n_cache.size(),0.0);
201 }
202
204
208 template<typename EG, typename LFSVC>
209 void onUnbindLFSV(const EG & eg, const LFSVC & lfsv_cache)
210 {
211 //global_rl_view.add(rl);
212 global_rl_view.unbind();
213 global_rl_view.commit();
214 }
215
216 template<typename IG, typename LFSVC>
217 void onUnbindLFSVInside(const IG & ig, const LFSVC & lfsv_cache)
218 {
219 //global_rl_view.add(rl);
220 global_rl_view.commit();
221 global_rl_view.unnbind();
222 }
223
224 template<typename IG, typename LFSVC>
225 void onUnbindLFSVOutside(const IG & ig,
226 const LFSVC & lfsv_s_cache,
227 const LFSVC & lfsv_n_cache)
228 {
229 //global_rn_view.add(rn);
230 global_rn_view.commit();
231 global_rn_view.unbind();
232 }
234
237 template<typename LFSUC>
238 void loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)
239 {
240 //global_sl_view.read(xl);
241 }
242 template<typename LFSUC>
243 void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
244 {
245 //global_sn_view.read(xn);
246 }
247 template<typename LFSUC>
248 void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
249 {
250 DUNE_THROW(Dune::NotImplemented,"No coupling lfsu available for ");
251 }
253
256
257 void postAssembly(const GFSU& gfsu, const GFSV& gfsv)
258 {
259 if(local_assembler.doPostProcessing())
260 Dune::PDELab::constrain_residual(local_assembler.testConstraints(),
261 global_rl_view.container());
262
263 }
264
266
269
270
277 template <typename EG> bool skipEntity(const EG &eg)
278 {
279 return localAssembler().skipEntity(eg) ||
280 (LocalAssembler::isNonOverlapping &&
281 eg.entity().partitionType() != Dune::InteriorEntity);
282 }
283
290 template<typename IG>
291 bool skipIntersection(const IG & ig)
292 {
293 return localAssembler().skipIntersection(ig);
294 }
295
296 template<typename EG, typename LFSUC, typename LFSVC>
297 void assembleUVVolume(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
298 {
299 //rl_view.setWeight(local_assembler.weight());
300 global_rl_view.setWeight(local_assembler.weight());
301 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolume>::
302 alpha_volume(lop,eg,lfsu_cache.localFunctionSpace(),global_sl_view,lfsv_cache.localFunctionSpace(),global_rl_view);
303 }
304
305 template<typename EG, typename LFSVC>
306 void assembleVVolume(const EG & eg, const LFSVC & lfsv_cache)
307 {
308 //rl_view.setWeight(local_assembler.weight());
309 global_rl_view.setWeight(local_assembler.weight());
310 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaVolume>::
311 lambda_volume(lop,eg,lfsv_cache.localFunctionSpace(),global_rl_view);
312 }
313
314 template<typename IG, typename LFSUC, typename LFSVC>
315 void assembleUVSkeleton(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
316 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
317 {
318 //rl_view.setWeight(local_assembler.weight());
319 //rn_view.setWeight(local_assembler.weight());
320 global_rl_view.setWeight(local_assembler.weight());
321 global_rn_view.setWeight(local_assembler.weight());
322 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaSkeleton>::
323 alpha_skeleton(lop,ig,
324 lfsu_s_cache.localFunctionSpace(),global_sl_view,lfsv_s_cache.localFunctionSpace(),
325 lfsu_n_cache.localFunctionSpace(),global_sn_view,lfsv_n_cache.localFunctionSpace(),
326 global_rl_view,global_rn_view);
327 }
328
329 template<typename IG, typename LFSVC>
330 void assembleVSkeleton(const IG & ig, const LFSVC & lfsv_s_cache, const LFSVC & lfsv_n_cache)
331 {
332 //rl_view.setWeight(local_assembler.weight());
333 //rn_view.setWeight(local_assembler.weight());
334 global_rl_view.setWeight(local_assembler.weight());
335 global_rn_view.setWeight(local_assembler.weight());
336 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaSkeleton>::
337 lambda_skeleton(lop, ig, lfsv_s_cache.localFunctionSpace(), lfsv_n_cache.localFunctionSpace(), global_rl_view, global_rn_view);
338 }
339
340 template<typename IG, typename LFSUC, typename LFSVC>
341 void assembleUVBoundary(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache)
342 {
343 //rl_view.setWeight(local_assembler.weight());
344 global_rl_view.setWeight(local_assembler.weight());
345 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaBoundary>::
346 alpha_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),global_sl_view,lfsv_s_cache.localFunctionSpace(),global_rl_view);
347 }
348
349 template<typename IG, typename LFSVC>
350 void assembleVBoundary(const IG & ig, const LFSVC & lfsv_s_cache)
351 {
352 //rl_view.setWeight(local_assembler.weight());
353 global_rl_view.setWeight(local_assembler.weight());
354 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaBoundary>::
355 lambda_boundary(lop,ig,lfsv_s_cache.localFunctionSpace(),global_rl_view);
356 }
357
358 template<typename IG, typename LFSUC, typename LFSVC>
359 static void assembleUVEnrichedCoupling(const IG & ig,
360 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
361 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache,
362 const LFSUC & lfsu_coupling_cache, const LFSVC & lfsv_coupling_cache)
363 {
364 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
365 }
366
367 template<typename IG, typename LFSVC>
368 static void assembleVEnrichedCoupling(const IG & ig,
369 const LFSVC & lfsv_s_cache,
370 const LFSVC & lfsv_n_cache,
371 const LFSVC & lfsv_coupling_cache)
372 {
373 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
374 }
375
376 template<typename EG, typename LFSUC, typename LFSVC>
377 void assembleUVVolumePostSkeleton(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
378 {
379 //rl_view.setWeight(local_assembler.weight());
380 global_rl_view.setWeight(local_assembler.weight());
381 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolumePostSkeleton>::
382 alpha_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),global_sl_view,lfsv_cache.localFunctionSpace(),global_rl_view);
383 }
384
385 template<typename EG, typename LFSVC>
386 void assembleVVolumePostSkeleton(const EG & eg, const LFSVC & lfsv_cache)
387 {
388 //rl_view.setWeight(local_assembler.weight());
389 global_rl_view.setWeight(local_assembler.weight());
390 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaVolumePostSkeleton>::
391 lambda_volume_post_skeleton(lop,eg,lfsv_cache.localFunctionSpace(),global_rl_view);
392 }
393
395
396 private:
399 const LocalAssembler & local_assembler;
400
402 const LOP & lop;
403
405 ResidualView global_rl_view;
406 ResidualView global_rn_view;
407
409 SolutionView global_sl_view;
410 SolutionView global_sn_view;
411
414 typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
415 typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
416
419
421 //typename ResidualVector::WeightedAccumulationView rl_view;
423 //typename ResidualVector::WeightedAccumulationView rn_view;
425
426 }; // End of class FastDGLocalResidualAssemblerEngine
427
428 }
429}
430#endif
Default exception for dummy implementations.
Definition: exceptions.hh:263
The fast DG local assembler engine for DUNE grids which assembles the residual vector.
Definition: residualengine.hh:24
LA::Traits::Residual Residual
The type of the residual vector.
Definition: residualengine.hh:40
LA::Traits::Solution Solution
The type of the solution vector.
Definition: residualengine.hh:44
void setSolution(const Solution &solution_)
Definition: residualengine.hh:148
LA LocalAssembler
The type of the wrapping local assembler.
Definition: residualengine.hh:34
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: residualengine.hh:158
bool requireSkeleton() const
Definition: residualengine.hh:98
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: residualengine.hh:133
void setResidual(Residual &residual_)
Definition: residualengine.hh:140
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: residualengine.hh:209
FastDGLocalResidualAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: residualengine.hh:64
LA::LFSU LFSU
The local function spaces.
Definition: residualengine.hh:48
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: residualengine.hh:238
bool skipIntersection(const IG &ig)
Definition: residualengine.hh:291
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: residualengine.hh:257
FastDGLocalResidualAssemblerEngine(const FastDGLocalResidualAssemblerEngine &other)
copy contructor
Definition: residualengine.hh:86
LA::LocalOperator LOP
The type of the local operator.
Definition: residualengine.hh:37
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: residualengine.hh:127
bool skipEntity(const EG &eg)
Definition: residualengine.hh:277
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: residualengine.hh:121
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:22
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
@ InteriorEntity
all interior entities
Definition: gridenums.hh:31
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:904
Dune namespace.
Definition: alignedallocator.hh:13
Definition: localfunctionspacetags.hh:54
Definition: localfunctionspacetags.hh:48
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)