1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_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>
26 template<
typename TrialConstra
intsContainer,
typename TestConstra
intsContainer>
27 bool needsConstraintsCaching(
const TrialConstraintsContainer& cu,
const TestConstraintsContainer& cv)
const
36 typedef typename LA::LocalOperator
LOP;
39 typedef typename LA::Traits::Residual
Residual;
40 typedef typename Residual::ElementType ResidualElement;
43 typedef typename LA::Traits::Solution
Solution;
44 typedef typename Solution::ElementType SolutionElement;
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;
54 typedef typename Solution::template ConstLocalView<LFSUCache> SolutionView;
55 typedef typename Residual::template LocalView<LFSVCache> ResidualView;
64 : local_assembler(local_assembler_),
65 lop(local_assembler_.localOperator()),
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(); }
97 return local_assembler;
101 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints&
trialConstraints()
const
107 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints&
testConstraints()
const
116 global_rl_view.attach(residual_);
117 global_rn_view.attach(residual_);
124 global_sl_view.attach(solution_);
125 global_sn_view.attach(solution_);
131 template<
typename EG,
typename LFSUC,
typename LFSVC>
132 void onBindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
134 global_sl_view.bind(lfsu_cache);
135 xl.
resize(lfsu_cache.size());
138 template<
typename EG,
typename LFSVC>
139 void onBindLFSV(
const EG & eg,
const LFSVC & lfsv_cache)
141 global_rl_view.bind(lfsv_cache);
142 rl.
assign(lfsv_cache.size(),0.0);
145 template<
typename IG,
typename LFSUC,
typename LFSVC>
146 void onBindLFSUVInside(
const IG & ig,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
148 global_sl_view.bind(lfsu_cache);
149 xl.
resize(lfsu_cache.size());
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)
157 global_sn_view.bind(lfsu_n_cache);
158 xn.
resize(lfsu_n_cache.size());
161 template<
typename IG,
typename LFSVC>
162 void onBindLFSVInside(
const IG & ig,
const LFSVC & lfsv_cache)
164 global_rl_view.bind(lfsv_cache);
165 rl.
assign(lfsv_cache.size(),0.0);
168 template<
typename IG,
typename LFSVC>
169 void onBindLFSVOutside(
const IG & ig,
170 const LFSVC & lfsv_s_cache,
171 const LFSVC & lfsv_n_cache)
173 global_rn_view.bind(lfsv_n_cache);
174 rn.
assign(lfsv_n_cache.size(),0.0);
182 template<
typename EG,
typename LFSVC>
185 global_rl_view.add(rl);
186 global_rl_view.commit();
189 template<
typename IG,
typename LFSVC>
190 void onUnbindLFSVInside(
const IG & ig,
const LFSVC & lfsv_cache)
192 global_rl_view.add(rl);
193 global_rl_view.commit();
196 template<
typename IG,
typename LFSVC>
197 void onUnbindLFSVOutside(
const IG & ig,
198 const LFSVC & lfsv_s_cache,
199 const LFSVC & lfsv_n_cache)
201 global_rn_view.add(rn);
202 global_rn_view.commit();
208 template<
typename LFSUC>
211 global_sl_view.read(xl);
213 template<
typename LFSUC>
214 void loadCoefficientsLFSUOutside(
const LFSUC & lfsu_n_cache)
216 global_sn_view.read(xn);
218 template<
typename LFSUC>
219 void loadCoefficientsLFSUCoupling(
const LFSUC & lfsu_c_cache)
230 if(local_assembler.doPostProcessing())
246 template<
typename EG>
258 template<
typename IG>
264 template<
typename EG,
typename LFSUC,
typename LFSVC>
265 void assembleUVVolume(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
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);
272 template<
typename EG,
typename LFSVC>
273 void assembleVVolume(
const EG & eg,
const LFSVC & lfsv_cache)
275 rl_view.
setWeight(local_assembler.weight());
276 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doLambdaVolume>::
277 lambda_volume(lop,eg,lfsv_cache.localFunctionSpace(),rl_view);
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)
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(),
293 template<
typename IG,
typename LFSVC>
294 void assembleVSkeleton(
const IG & ig,
const LFSVC & lfsv_s_cache,
const LFSVC & lfsv_n_cache)
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);
302 template<
typename IG,
typename LFSUC,
typename LFSVC>
303 void assembleUVBoundary(
const IG & ig,
const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache)
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);
310 template<
typename IG,
typename LFSVC>
311 void assembleVBoundary(
const IG & ig,
const LFSVC & lfsv_s_cache)
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);
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)
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)
336 template<
typename EG,
typename LFSUC,
typename LFSVC>
337 void assembleUVVolumePostSkeleton(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
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);
344 template<
typename EG,
typename LFSVC>
345 void assembleVVolumePostSkeleton(
const EG & eg,
const LFSVC & lfsv_cache)
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);
363 ResidualView global_rl_view;
364 ResidualView global_rn_view;
367 SolutionView global_sl_view;
368 SolutionView global_sn_view;
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