1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH 2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH 26 template<
typename TrialConstra
intsContainer,
typename TestConstra
intsContainer>
36 typedef typename LA::LocalOperator
LOP;
39 typedef typename LA::Traits::Residual
Residual;
43 typedef typename LA::Traits::Solution
Solution;
47 typedef typename LA::LFSU
LFSU;
49 typedef typename LFSU::Traits::GridFunctionSpace
GFSU;
50 typedef typename LA::LFSV
LFSV;
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_), lop(local_assembler_.lop),
72 {
return ( local_assembler.doAlphaSkeleton() || local_assembler.doLambdaSkeleton() ); }
74 {
return local_assembler.doSkeletonTwoSided(); }
76 {
return local_assembler.doAlphaVolume(); }
78 {
return local_assembler.doLambdaVolume(); }
80 {
return local_assembler.doAlphaSkeleton(); }
82 {
return local_assembler.doLambdaSkeleton(); }
84 {
return local_assembler.doAlphaBoundary(); }
86 {
return local_assembler.doLambdaBoundary(); }
88 {
return local_assembler.doAlphaVolumePostSkeleton(); }
90 {
return local_assembler.doLambdaVolumePostSkeleton(); }
97 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints&
trialConstraints()
const 103 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints&
testConstraints()
const 111 global_rl_view.attach(residual_);
112 global_rn_view.attach(residual_);
118 global_sl_view.attach(solution_);
119 global_sn_view.attach(solution_);
125 template<
typename EG,
typename LFSUC,
typename LFSVC>
126 void onBindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache){
127 global_sl_view.bind(lfsu_cache);
128 xl.
resize(lfsu_cache.size());
131 template<
typename EG,
typename LFSVC>
133 global_rl_view.bind(lfsv_cache);
134 rl.
assign(lfsv_cache.size(),0.0);
137 template<
typename IG,
typename LFSUC,
typename LFSVC>
139 global_sl_view.bind(lfsu_cache);
140 xl.
resize(lfsu_cache.size());
143 template<
typename IG,
typename LFSUC,
typename LFSVC>
145 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
146 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
148 global_sn_view.bind(lfsu_n_cache);
149 xn.
resize(lfsu_n_cache.size());
152 template<
typename IG,
typename LFSVC>
154 global_rl_view.bind(lfsv_cache);
155 rl.
assign(lfsv_cache.size(),0.0);
158 template<
typename IG,
typename LFSVC>
160 const LFSVC & lfsv_s_cache,
161 const LFSVC & lfsv_n_cache)
163 global_rn_view.bind(lfsv_n_cache);
164 rn.
assign(lfsv_n_cache.size(),0.0);
172 template<
typename EG,
typename LFSVC>
174 global_rl_view.add(rl);
175 global_rl_view.commit();
178 template<
typename IG,
typename LFSVC>
180 global_rl_view.add(rl);
181 global_rl_view.commit();
184 template<
typename IG,
typename LFSVC>
186 const LFSVC & lfsv_s_cache,
187 const LFSVC & lfsv_n_cache)
189 global_rn_view.add(rn);
190 global_rn_view.commit();
196 template<
typename LFSUC>
198 global_sl_view.read(xl);
200 template<
typename LFSUC>
202 global_sn_view.read(xn);
204 template<
typename LFSUC>
206 {DUNE_THROW(Dune::NotImplemented,
"No coupling lfsu available for ");}
213 if(local_assembler.doPostProcessing)
230 template<
typename EG>
233 return LocalAssembler::isNonOverlapping && eg.entity().partitionType() != Dune::InteriorEntity;
236 template<
typename EG,
typename LFSUC,
typename LFSVC>
239 rl_view.
setWeight(local_assembler.weight);
241 alpha_volume(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),rl_view);
244 template<
typename EG,
typename LFSVC>
247 rl_view.
setWeight(local_assembler.weight);
249 lambda_volume(lop,eg,lfsv_cache.localFunctionSpace(),rl_view);
252 template<
typename IG,
typename LFSUC,
typename LFSVC>
254 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
256 rl_view.
setWeight(local_assembler.weight);
257 rn_view.
setWeight(local_assembler.weight);
260 lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),
261 lfsu_n_cache.localFunctionSpace(),xn,lfsv_n_cache.localFunctionSpace(),
265 template<
typename IG,
typename LFSVC>
268 rl_view.
setWeight(local_assembler.weight);
269 rn_view.
setWeight(local_assembler.weight);
271 lambda_skeleton(lop, ig, lfsv_s_cache.localFunctionSpace(), lfsv_n_cache.localFunctionSpace(), rl_view, rn_view);
274 template<
typename IG,
typename LFSUC,
typename LFSVC>
277 rl_view.
setWeight(local_assembler.weight);
279 alpha_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),rl_view);
282 template<
typename IG,
typename LFSVC>
285 rl_view.
setWeight(local_assembler.weight);
290 template<
typename IG,
typename LFSUC,
typename LFSVC>
292 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
293 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache,
294 const LFSUC & lfsu_coupling_cache,
const LFSVC & lfsv_coupling_cache)
295 {DUNE_THROW(Dune::NotImplemented,
"Assembling of coupling spaces is not implemented for ");}
297 template<
typename IG,
typename LFSVC>
299 const LFSVC & lfsv_s_cache,
300 const LFSVC & lfsv_n_cache,
301 const LFSVC & lfsv_coupling_cache)
302 {DUNE_THROW(Dune::NotImplemented,
"Assembling of coupling spaces is not implemented for ");}
304 template<
typename EG,
typename LFSUC,
typename LFSVC>
307 rl_view.
setWeight(local_assembler.weight);
312 template<
typename EG,
typename LFSVC>
315 rl_view.
setWeight(local_assembler.weight);
325 const LocalAssembler & local_assembler;
331 ResidualView global_rl_view;
332 ResidualView global_rn_view;
335 SolutionView global_sl_view;
336 SolutionView global_sn_view;
364 #endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_RESIDUALENGINE_HH static void alpha_volume_post_skeleton(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r)
Definition: callswitch.hh:95
const IG & ig
Definition: constraints.hh:148
bool needsConstraintsCaching(const TrialConstraintsContainer &cu, const TestConstraintsContainer &cv) const
Definition: default/residualengine.hh:27
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: default/residualengine.hh:103
DefaultLocalResidualAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: default/residualengine.hh:63
bool requireVBoundary() const
Definition: default/residualengine.hh:85
void assembleVVolumePostSkeleton(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:313
LA::Traits::Solution Solution
The type of the solution vector.
Definition: default/residualengine.hh:43
bool requireUVBoundary() const
Definition: default/residualengine.hh:83
static void lambda_skeleton(const LA &la, const IG &ig, const LFSV &lfsv_s, const LFSV &lfsv_n, R &r_s, R &r_n)
Definition: callswitch.hh:121
void setSolution(const Solution &solution_)
Definition: default/residualengine.hh:117
void assembleUVVolumePostSkeleton(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:305
bool requireVSkeleton() const
Definition: default/residualengine.hh:81
static void 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)
Definition: default/residualengine.hh:291
void onUnbindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: default/residualengine.hh:185
LA::LocalOperator LOP
The type of the local operator.
Definition: default/residualengine.hh:36
LA::LFSU LFSU
The local function spaces.
Definition: default/residualengine.hh:47
Definition: localfunctionspacetags.hh:54
static void alpha_volume(const LA &la, const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, R &r)
Definition: callswitch.hh:91
bool requireSkeletonTwoSided() const
Definition: default/residualengine.hh:73
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: default/residualengine.hh:97
void assembleUVVolume(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:237
Solution::ElementType SolutionElement
Definition: default/residualengine.hh:44
static void lambda_volume(const LA &la, const EG &eg, const LFSV &lfsv, R &r)
Definition: callswitch.hh:113
void assign(size_type size, const T &value)
Resize the container to size and assign the passed value to all entries.
Definition: localvector.hh:257
void onUnbindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:173
void setResidual(Residual &residual_)
Definition: default/residualengine.hh:110
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:126
Definition: localfunctionspacetags.hh:48
static void alpha_skeleton(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, const LFSU &lfsu_n, const X &x_n, const LFSV &lfsv_n, R &r_s, R &r_n)
Definition: callswitch.hh:99
LA::LFSV LFSV
Definition: default/residualengine.hh:50
Solution::template ConstLocalView< LFSUCache > SolutionView
Definition: default/residualengine.hh:54
void setWeight(weight_type weight)
Resets the weighting coefficient of the view.
Definition: localvector.hh:72
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: default/residualengine.hh:94
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: default/residualengine.hh:212
static void lambda_boundary(const LA &la, const IG &ig, const LFSV &lfsv, R &r)
Definition: callswitch.hh:127
LFSV::Traits::GridFunctionSpace GFSV
Definition: default/residualengine.hh:52
void onBindLFSUVInside(const IG &ig, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:138
static void assembleVEnrichedCoupling(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache, const LFSVC &lfsv_coupling_cache)
Definition: default/residualengine.hh:298
void onUnbindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:179
bool requireUVVolumePostSkeleton() const
Definition: default/residualengine.hh:87
bool requireVVolume() const
Definition: default/residualengine.hh:77
LFSU::Traits::GridFunctionSpace GFSU
Definition: default/residualengine.hh:49
void onBindLFSUVOutside(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache, const LFSUC &lfsu_n_cache, const LFSVC &lfsv_n_cache)
Definition: default/residualengine.hh:144
void assembleVBoundary(const IG &ig, const LFSVC &lfsv_s_cache)
Definition: default/residualengine.hh:283
bool requireUVSkeleton() const
Definition: default/residualengine.hh:79
void loadCoefficientsLFSUInside(const LFSUC &lfsu_s_cache)
Definition: default/residualengine.hh:197
Residual::template LocalView< LFSVCache > ResidualView
Definition: default/residualengine.hh:55
void onBindLFSVOutside(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: default/residualengine.hh:159
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
bool requireUVVolume() const
Definition: default/residualengine.hh:75
LA LocalAssembler
The type of the wrapping local assembler.
Definition: default/residualengine.hh:33
static void alpha_boundary(const LA &la, const IG &ig, const LFSU &lfsu_s, const X &x_s, const LFSV &lfsv_s, R &r_s)
Definition: callswitch.hh:106
bool requireSkeleton() const
Definition: default/residualengine.hh:71
LA::LFSVCache LFSVCache
Definition: default/residualengine.hh:51
WeightedVectorAccumulationView< LocalVector > WeightedAccumulationView
An accumulate-only view of this container that automatically applies a weight to all contributions...
Definition: localvector.hh:198
The local assembler engine for DUNE grids which assembles the residual vector.
Definition: default/residualengine.hh:21
static void lambda_volume_post_skeleton(const LA &la, const EG &eg, const LFSV &lfsv, R &r)
Definition: callswitch.hh:117
bool assembleCell(const EG &eg)
Definition: default/residualengine.hh:231
void assembleVVolume(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:245
void assembleUVBoundary(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache)
Definition: default/residualengine.hh:275
void constrain_residual(const CG &cg, XG &xg)
transform residual into transformed basis: r -> r~
Definition: constraints.hh:906
void onBindLFSVInside(const IG &ig, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:153
void assembleUVSkeleton(const IG &ig, const LFSUC &lfsu_s_cache, const LFSVC &lfsv_s_cache, const LFSUC &lfsu_n_cache, const LFSVC &lfsv_n_cache)
Definition: default/residualengine.hh:253
void resize(size_type size)
Resize the container.
Definition: localvector.hh:251
bool requireVVolumePostSkeleton() const
Definition: default/residualengine.hh:89
void assembleVSkeleton(const IG &ig, const LFSVC &lfsv_s_cache, const LFSVC &lfsv_n_cache)
Definition: default/residualengine.hh:266
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:21
LA::Traits::Residual Residual
The type of the residual vector.
Definition: default/residualengine.hh:39
Residual::ElementType ResidualElement
Definition: default/residualengine.hh:40
void loadCoefficientsLFSUCoupling(const LFSUC &lfsu_c_cache)
Definition: default/residualengine.hh:205
void loadCoefficientsLFSUOutside(const LFSUC &lfsu_n_cache)
Definition: default/residualengine.hh:201
void onBindLFSV(const EG &eg, const LFSVC &lfsv_cache)
Definition: default/residualengine.hh:132
LA::LFSUCache LFSUCache
Definition: default/residualengine.hh:48