DUNE PDELab (2.8)

jacobianengine.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_FASTDG_JACOBIANENGINE_HH
2#define DUNE_PDELAB_GRIDOPERATOR_FASTDG_JACOBIANENGINE_HH
3
4#include <dune/pdelab/constraints/common/constraints.hh>
6#include <dune/pdelab/gridoperator/common/localmatrix.hh>
7#include <dune/pdelab/gridoperator/common/diagonallocalmatrix.hh>
8#include <dune/pdelab/gridoperator/common/gridoperatorutilities.hh>
9#include <dune/pdelab/gridoperator/common/assemblerutilities.hh>
10#include <dune/pdelab/gridoperator/common/localassemblerenginebase.hh>
11#include <dune/pdelab/localoperator/callswitch.hh>
12#include <dune/pdelab/localoperator/flags.hh>
13
14namespace Dune{
15 namespace PDELab{
16
24 template<typename LA>
27 {
28 public:
29
30 template<typename TrialConstraintsContainer, typename TestConstraintsContainer>
31 bool needsConstraintsCaching(const TrialConstraintsContainer& cu, const TestConstraintsContainer& cv)
32 {
33 return cu.containsNonDirichletConstraints() || cv.containsNonDirichletConstraints();
34 }
35
37 typedef LA LocalAssembler;
38
40 typedef typename LA::LocalOperator LOP;
41
43 typedef typename LA::LFSU LFSU;
44 typedef typename LA::LFSUCache LFSUCache;
45 typedef typename LFSU::Traits::GridFunctionSpace GFSU;
46 typedef typename LA::LFSV LFSV;
47 typedef typename LA::LFSVCache LFSVCache;
48 typedef typename LFSV::Traits::GridFunctionSpace GFSV;
49
51 typedef typename LA::Traits::Jacobian Jacobian;
52 typedef typename Jacobian::ElementType JacobianElement;
53 typedef typename Jacobian::template AliasedLocalView<LFSVCache,LFSUCache> JacobianView;
54
56 typedef typename LA::Traits::Solution Solution;
57 typedef typename Solution::ElementType SolutionElement;
58 typedef typename Solution::template ConstAliasedLocalView<LFSUCache> SolutionView;
59
67 : local_assembler(local_assembler_),
68 lop(local_assembler_.localOperator())
69 {}
70
72
87 local_assembler(other.local_assembler), lop(other.lop),
88 global_s_s_view(other.global_s_s_view),
89 global_s_n_view(other.global_s_n_view),
90 global_a_ss_view(other.global_a_ss_view),
91 global_a_sn_view(other.global_a_sn_view),
92 global_a_ns_view(other.global_a_ns_view),
93 global_a_nn_view(other.global_a_nn_view)
94 { }
95
98 bool requireSkeleton() const
99 { return local_assembler.doAlphaSkeleton(); }
100 bool requireSkeletonTwoSided() const
101 { return local_assembler.doSkeletonTwoSided(); }
102 bool requireUVVolume() const
103 { return local_assembler.doAlphaVolume(); }
104 bool requireUVSkeleton() const
105 { return local_assembler.doAlphaSkeleton(); }
106 bool requireUVBoundary() const
107 { return local_assembler.doAlphaBoundary(); }
108 bool requireUVVolumePostSkeleton() const
109 { return local_assembler.doAlphaVolumePostSkeleton(); }
111
114 {
115 return local_assembler;
116 }
117
119 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints& trialConstraints() const
120 {
121 return localAssembler().trialConstraints();
122 }
123
125 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints& testConstraints() const
126 {
127 return localAssembler().testConstraints();
128 }
129
132 void setJacobian(Jacobian & jacobian_)
133 {
134 global_a_ss_view.attach(jacobian_);
135 global_a_sn_view.attach(jacobian_);
136 global_a_ns_view.attach(jacobian_);
137 global_a_nn_view.attach(jacobian_);
138 }
139
142 void setSolution(const Solution & solution_)
143 {
144 global_s_s_view.attach(solution_);
145 global_s_n_view.attach(solution_);
146 }
147
151 template<typename EG, typename LFSUC, typename LFSVC>
152 void onBindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
153 {
154 global_s_s_view.bind(lfsu_cache);
155 global_a_ss_view.bind(lfsv_cache,lfsu_cache);
156 }
157
158 template<typename IG, typename LFSUC, typename LFSVC>
159 void onBindLFSUVOutside(const IG & ig,
160 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
161 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
162 {
163 global_s_n_view.bind(lfsu_n_cache);
164 global_a_sn_view.bind(lfsv_s_cache,lfsu_n_cache);
165 global_a_ns_view.bind(lfsv_n_cache,lfsu_s_cache);
166 global_a_nn_view.bind(lfsv_n_cache,lfsu_n_cache);
167 }
168
170
174 template<typename EG, typename LFSUC, typename LFSVC>
175 void onUnbindLFSUV(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
176 {
177 global_a_ss_view.unbind();
178 }
179
180 template<typename IG, typename LFSUC, typename LFSVC>
181 void onUnbindLFSUVOutside(const IG & ig,
182 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
183 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
184 {
185 global_a_sn_view.unbind();
186 global_a_ns_view.unbind();
187 global_a_nn_view.unbind();
188 }
189
191
194 template<typename LFSUC>
195 void loadCoefficientsLFSUInside(const LFSUC & lfsu_cache)
196 {}
197 template<typename LFSUC>
198 void loadCoefficientsLFSUOutside(const LFSUC & lfsu_n_cache)
199 {}
200 template<typename LFSUC>
201 void loadCoefficientsLFSUCoupling(const LFSUC & lfsu_c_cache)
202 {
203 DUNE_THROW(Dune::NotImplemented,"No coupling lfsu_cache available for ");
204 }
206
209 void postAssembly(const GFSU& gfsu, const GFSV& gfsv)
210 {
211 Jacobian& jacobian = global_a_ss_view.container();
212 global_s_s_view.detach();
213 global_s_n_view.detach();
214 global_a_ss_view.detach();
215 global_a_sn_view.detach();
216 global_a_ns_view.detach();
217 global_a_nn_view.detach();
218
219 if(local_assembler.doPostProcessing())
220 local_assembler.handle_dirichlet_constraints(gfsv,jacobian);
221 }
223
226
233 template <typename EG> bool skipEntity(const EG &eg)
234 {
235 return localAssembler().skipEntity(eg) ||
236 (LocalAssembler::isNonOverlapping &&
237 eg.entity().partitionType() != Dune::InteriorEntity);
238 }
239
246 template<typename IG>
247 bool skipIntersection(const IG & ig)
248 {
249 return localAssembler().skipIntersection(ig);
250 }
251
252 template<typename EG, typename LFSUC, typename LFSVC>
253 void assembleUVVolume(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
254 {
255 global_a_ss_view.setWeight(local_assembler.weight());
256 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolume>::
257 jacobian_volume(lop,eg,lfsu_cache.localFunctionSpace(),global_s_s_view,lfsv_cache.localFunctionSpace(),global_a_ss_view);
258 }
259
260 template<typename IG, typename LFSUC, typename LFSVC>
261 void assembleUVSkeleton(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
262 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache)
263 {
264 global_a_ss_view.setWeight(local_assembler.weight());
265 global_a_sn_view.setWeight(local_assembler.weight());
266 global_a_ns_view.setWeight(local_assembler.weight());
267 global_a_nn_view.setWeight(local_assembler.weight());
268 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaSkeleton>::
269 jacobian_skeleton(lop,ig,
270 lfsu_s_cache.localFunctionSpace(),global_s_s_view,lfsv_s_cache.localFunctionSpace(),
271 lfsu_n_cache.localFunctionSpace(),global_s_n_view,lfsv_n_cache.localFunctionSpace(),
272 global_a_ss_view, global_a_sn_view,
273 global_a_ns_view, global_a_nn_view);
274 }
275
276 template<typename IG, typename LFSUC, typename LFSVC>
277 void assembleUVBoundary(const IG & ig, const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache)
278 {
279 global_a_ss_view.setWeight(local_assembler.weight());
280 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaBoundary>::
281 jacobian_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),global_s_s_view,lfsv_s_cache.localFunctionSpace(),global_a_ss_view);
282 }
283
284 template<typename IG, typename LFSUC, typename LFSVC>
285 static void assembleUVEnrichedCoupling(const IG & ig,
286 const LFSUC & lfsu_s_cache, const LFSVC & lfsv_s_cache,
287 const LFSUC & lfsu_n_cache, const LFSVC & lfsv_n_cache,
288 const LFSUC & lfsu_coupling_cache, const LFSVC & lfsv_coupling_cache)
289 {
290 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
291 }
292
293 template<typename IG, typename LFSVC>
294 static void assembleVEnrichedCoupling(const IG & ig,
295 const LFSVC & lfsv_s_cache,
296 const LFSVC & lfsv_n_cache,
297 const LFSVC & lfsv_coupling_cache)
298 {
299 DUNE_THROW(Dune::NotImplemented,"Assembling of coupling spaces is not implemented for ");
300 }
301
302 template<typename EG, typename LFSUC, typename LFSVC>
303 void assembleUVVolumePostSkeleton(const EG & eg, const LFSUC & lfsu_cache, const LFSVC & lfsv_cache)
304 {
305 global_a_ss_view.setWeight(local_assembler.weight());
306 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolumePostSkeleton>::
307 jacobian_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),global_s_s_view,lfsv_cache.localFunctionSpace(),global_a_ss_view);
308 }
309
311
312 private:
315 const LocalAssembler & local_assembler;
316
318 const LOP & lop;
319
321 SolutionView global_s_s_view;
322 SolutionView global_s_n_view;
323
325 JacobianView global_a_ss_view;
326 JacobianView global_a_sn_view;
327 JacobianView global_a_ns_view;
328 JacobianView global_a_nn_view;
329
332 typedef Dune::PDELab::TrialSpaceTag LocalTrialSpaceTag;
333 typedef Dune::PDELab::TestSpaceTag LocalTestSpaceTag;
334
336 typedef typename std::conditional<
337 std::is_base_of<
338 lop::DiagonalJacobian,
339 LOP
340 >::value,
343 >::type JacobianMatrix;
344
346
347 }; // End of class FastDGLocalJacobianAssemblerEngine
348
349 }
350}
351#endif
Default exception for dummy implementations.
Definition: exceptions.hh:261
A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpac...
Definition: diagonallocalmatrix.hh:29
The fast DG local assembler engine for DUNE grids which assembles the jacobian matrix.
Definition: jacobianengine.hh:27
bool requireSkeleton() const
Definition: jacobianengine.hh:98
LA::Traits::Solution Solution
The type of the solution vector.
Definition: jacobianengine.hh:56
bool skipEntity(const EG &eg)
Definition: jacobianengine.hh:233
void setSolution(const Solution &solution_)
Definition: jacobianengine.hh:142
FastDGLocalJacobianAssemblerEngine(const FastDGLocalJacobianAssemblerEngine &other)
copy contructor
Definition: jacobianengine.hh:86
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: jacobianengine.hh:113
LA::LFSU LFSU
The local function spaces.
Definition: jacobianengine.hh:43
void setJacobian(Jacobian &jacobian_)
Definition: jacobianengine.hh:132
LA LocalAssembler
The type of the wrapping local assembler.
Definition: jacobianengine.hh:37
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: jacobianengine.hh:125
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: jacobianengine.hh:119
bool skipIntersection(const IG &ig)
Definition: jacobianengine.hh:247
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: jacobianengine.hh:209
LA::LocalOperator LOP
The type of the local operator.
Definition: jacobianengine.hh:40
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianengine.hh:152
void onUnbindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianengine.hh:175
FastDGLocalJacobianAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: jacobianengine.hh:66
void loadCoefficientsLFSUInside(const LFSUC &lfsu_cache)
Definition: jacobianengine.hh:195
LA::Traits::Jacobian Jacobian
The type of the jacobian matrix.
Definition: jacobianengine.hh:51
Base class for LocalAssemblerEngine implementations to avoid boilerplate code.
Definition: localassemblerenginebase.hh:22
A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpac...
Definition: localmatrix.hh:184
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
@ InteriorEntity
all interior entities
Definition: gridenums.hh:29
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 22, 23:30, 2024)