1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
3 
4 #include <dune/pdelab/gridfunctionspace/localvector.hh>
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 
10 namespace Dune{
11   namespace PDELab{
12 
13     /**
14        \brief The local assembler engine for DUNE grids which
15        assembles the residual vector
16 
17        \tparam LA The local assembler
18 
19     */
20     template<typename LA>
21     class DefaultLocalResidualAssemblerEngine
22       : public LocalAssemblerEngineBase
23     {
24     public:
25 
26       template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
needsConstraintsCaching(const TrialConstraintsContainer & cu,const TestConstraintsContainer & cv) const27       bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv) const
28       {
29         return false;
30       }
31 
32       //! The type of the wrapping local assembler
33       typedef LA LocalAssembler;
34 
35       //! The type of the local operator
36       typedef typename LA::LocalOperator LOP;
37 
38       //! The type of the residual vector
39       typedef typename LA::Traits::Residual Residual;
40       typedef typename Residual::ElementType ResidualElement;
41 
42       //! The type of the solution vector
43       typedef typename LA::Traits::Solution Solution;
44       typedef typename Solution::ElementType SolutionElement;
45 
46       //! The local function spaces
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 
57       /**
58          \brief Constructor
59 
60          \param [in] local_assembler_ The local assembler object which
61          creates this engine
62       */
DefaultLocalResidualAssemblerEngine(const LocalAssembler & local_assembler_)63       DefaultLocalResidualAssemblerEngine(const LocalAssembler & local_assembler_)
64         : local_assembler(local_assembler_),
65           lop(local_assembler_.localOperator()),
66           rl_view(rl,1.0),
67           rn_view(rn,1.0)
68       {}
69 
70       //! Query methods for the global grid assembler
71       //! @{
requireSkeleton() const72       bool requireSkeleton() const
73       { return ( local_assembler.doAlphaSkeleton() || local_assembler.doLambdaSkeleton() ); }
requireSkeletonTwoSided() const74       bool requireSkeletonTwoSided() const
75       { return local_assembler.doSkeletonTwoSided(); }
requireUVVolume() const76       bool requireUVVolume() const
77       { return local_assembler.doAlphaVolume(); }
requireVVolume() const78       bool requireVVolume() const
79       { return local_assembler.doLambdaVolume(); }
requireUVSkeleton() const80       bool requireUVSkeleton() const
81       { return local_assembler.doAlphaSkeleton(); }
requireVSkeleton() const82       bool requireVSkeleton() const
83       { return local_assembler.doLambdaSkeleton(); }
requireUVBoundary() const84       bool requireUVBoundary() const
85       { return local_assembler.doAlphaBoundary(); }
requireVBoundary() const86       bool requireVBoundary() const
87       { return local_assembler.doLambdaBoundary(); }
requireUVVolumePostSkeleton() const88       bool requireUVVolumePostSkeleton() const
89       { return local_assembler.doAlphaVolumePostSkeleton(); }
requireVVolumePostSkeleton() const90       bool requireVVolumePostSkeleton() const
91       { return local_assembler.doLambdaVolumePostSkeleton(); }
92       //! @}
93 
94       //! Public access to the wrapping local assembler
localAssembler() const95       const LocalAssembler & localAssembler() const
96       {
97         return local_assembler;
98       }
99 
100       //! Trial space constraints
trialConstraints() const101       const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
102       {
103         return localAssembler().trialConstraints();
104       }
105 
106       //! Test space constraints
testConstraints() const107       const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
108       {
109         return localAssembler().testConstraints();
110       }
111 
112       //! Set current residual vector. Should be called prior to
113       //! assembling.
setResidual(Residual & residual_)114       void setResidual(Residual & residual_)
115       {
116         global_rl_view.attach(residual_);
117         global_rn_view.attach(residual_);
118       }
119 
120       //! Set current solution vector. Should be called prior to
121       //! assembling.
setSolution(const Solution & solution_)122       void setSolution(const Solution & solution_)
123       {
124         global_sl_view.attach(solution_);
125         global_sn_view.attach(solution_);
126       }
127 
128       //! Called immediately after binding of local function space in
129       //! global assembler.
130       //! @{
131       template<typename EG, typename LFSUC, typename LFSVC>
onBindLFSUV(const EG & eg,const LFSUC & lfsu_cache,const LFSVC & lfsv_cache)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>
onBindLFSV(const EG & eg,const LFSVC & lfsv_cache)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>
onBindLFSUVInside(const IG & ig,const LFSUC & lfsu_cache,const LFSVC & lfsv_cache)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>
onBindLFSUVOutside(const IG & ig,const LFSUC & lfsu_s_cache,const LFSVC & lfsv_s_cache,const LFSUC & lfsu_n_cache,const LFSVC & lfsv_n_cache)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>
onBindLFSVInside(const IG & ig,const LFSVC & lfsv_cache)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>
onBindLFSVOutside(const IG & ig,const LFSVC & lfsv_s_cache,const LFSVC & lfsv_n_cache)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 
177       //! @}
178 
179       //! Called when the local function space is about to be rebound or
180       //! discarded
181       //! @{
182       template<typename EG, typename LFSVC>
onUnbindLFSV(const EG & eg,const LFSVC & lfsv_cache)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>
onUnbindLFSVInside(const IG & ig,const LFSVC & lfsv_cache)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>
onUnbindLFSVOutside(const IG & ig,const LFSVC & lfsv_s_cache,const LFSVC & lfsv_n_cache)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       }
204       //! @}
205 
206       //! Methods for loading of the local function's coefficients
207       //! @{
208       template<typename LFSUC>
loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)209       void loadCoefficientsLFSUInside(const LFSUC & lfsu_s_cache)
210       {
211         global_sl_view.read(xl);
212       }
213       template<typename LFSUC>
loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)214       void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
215       {
216         global_sn_view.read(xn);
217       }
218       template<typename LFSUC>
loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)219       void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
220       {
221         DUNE_THROW(Dune::NotImplemented,"No coupling lfsu available for ");
222       }
223       //! @}
224 
225       //! Notifier functions, called immediately before and after assembling
226       //! @{
227 
postAssembly(const GFSU & gfsu,const GFSV & gfsv)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 
235       //! @}
236 
237       //! Assembling methods
238       //! @{
239 
240       /** Assemble on a given cell without function spaces.
241 
242           \return If true, the assembling for this cell is assumed to
243           be complete and the assembler continues with the next grid
244           cell.
245        */
246       template<typename EG>
skipEntity(const EG & eg)247       bool skipEntity(const EG & eg)
248       {
249         return localAssembler().skipEntity(eg);
250       }
251 
252       /** Assemble on a given intersection without function spaces.
253 
254           \return If true, the assembling for this intersection is assumed to
255           be complete and the assembler continues with the next grid
256           intersection.
257        */
258       template<typename IG>
skipIntersection(const IG & ig)259       bool skipIntersection(const IG & ig)
260       {
261         return localAssembler().skipIntersection(ig);
262       }
263 
264       template<typename EG, typename LFSUC, typename LFSVC>
assembleUVVolume(const EG & eg,const LFSUC & lfsu_cache,const LFSVC & lfsv_cache)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>
assembleVVolume(const EG & eg,const LFSVC & lfsv_cache)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>
assembleUVSkeleton(const IG & ig,const LFSUC & lfsu_s_cache,const LFSVC & lfsv_s_cache,const LFSUC & lfsu_n_cache,const LFSVC & lfsv_n_cache)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>
assembleVSkeleton(const IG & ig,const LFSVC & lfsv_s_cache,const LFSVC & lfsv_n_cache)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>
assembleUVBoundary(const IG & ig,const LFSUC & lfsu_s_cache,const LFSVC & lfsv_s_cache)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>
assembleVBoundary(const IG & ig,const LFSVC & lfsv_s_cache)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>
assembleUVEnrichedCoupling(const IG & ig,const LFSUC & lfsu_s_cache,const LFSVC & lfsv_s_cache,const LFSUC & lfsu_n_cache,const LFSVC & lfsv_n_cache,const LFSUC & lfsu_coupling_cache,const LFSVC & lfsv_coupling_cache)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>
assembleVEnrichedCoupling(const IG & ig,const LFSVC & lfsv_s_cache,const LFSVC & lfsv_n_cache,const LFSVC & lfsv_coupling_cache)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>
assembleUVVolumePostSkeleton(const EG & eg,const LFSUC & lfsu_cache,const LFSVC & lfsv_cache)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>
assembleVVolumePostSkeleton(const EG & eg,const LFSVC & lfsv_cache)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 
352       //! @}
353 
354     private:
355       //! Reference to the wrapping local assembler object which
356       //! constructed this engine
357       const LocalAssembler & local_assembler;
358 
359       //! Reference to the local operator
360       const LOP & lop;
361 
362       //! Pointer to the current residual vector in which to assemble
363       ResidualView global_rl_view;
364       ResidualView global_rn_view;
365 
366       //! Pointer to the current residual vector in which to assemble
367       SolutionView global_sl_view;
368       SolutionView global_sn_view;
369 
370       //! The local vectors and matrices as required for assembling
371       //! @{
372       typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
373       typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
374 
375       typedef Dune::PDELab::LocalVector<SolutionElement, LocalTrialSpaceTag> SolutionVector;
376       typedef Dune::PDELab::LocalVector<ResidualElement, LocalTestSpaceTag> ResidualVector;
377 
378       //! Inside local coefficients
379       SolutionVector xl;
380       //! Outside local coefficients
381       SolutionVector xn;
382       //! Inside local residual
383       ResidualVector rl;
384       //! Outside local residual
385       ResidualVector rn;
386       //! Inside local residual weighted view
387       typename ResidualVector::WeightedAccumulationView rl_view;
388       //! Outside local residual weighted view
389       typename ResidualVector::WeightedAccumulationView rn_view;
390       //! @}
391 
392     }; // End of class DefaultLocalResidualAssemblerEngine
393 
394   }
395 }
396 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
397