1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_JACOBIANENGINE_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_JACOBIANENGINE_HH
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/assemblerutilities.hh>
9#include <dune/pdelab/gridoperator/common/localassemblerenginebase.hh>
10#include <dune/pdelab/localoperator/callswitch.hh>
11#include <dune/pdelab/localoperator/flags.hh>
29 template<
typename TrialConstra
intsContainer,
typename TestConstra
intsContainer>
30 bool needsConstraintsCaching(
const TrialConstraintsContainer& cu,
const TestConstraintsContainer& cv)
32 return cu.containsNonDirichletConstraints() || cv.containsNonDirichletConstraints();
39 typedef typename LA::LocalOperator
LOP;
42 typedef typename LA::LFSU
LFSU;
43 typedef typename LA::LFSUCache LFSUCache;
44 typedef typename LFSU::Traits::GridFunctionSpace GFSU;
45 typedef typename LA::LFSV LFSV;
46 typedef typename LA::LFSVCache LFSVCache;
47 typedef typename LFSV::Traits::GridFunctionSpace GFSV;
50 typedef typename LA::Traits::Jacobian
Jacobian;
51 typedef typename Jacobian::ElementType JacobianElement;
52 typedef typename Jacobian::template LocalView<LFSVCache,LFSUCache> JacobianView;
55 typedef typename LA::Traits::Solution
Solution;
56 typedef typename Solution::ElementType SolutionElement;
57 typedef typename Solution::template ConstLocalView<LFSUCache> SolutionView;
66 : local_assembler(local_assembler_),
67 lop(local_assembler_.localOperator()),
69 al_sn_view(al_sn,1.0),
70 al_ns_view(al_ns,1.0),
77 {
return local_assembler.doAlphaSkeleton(); }
78 bool requireSkeletonTwoSided()
const
79 {
return local_assembler.doSkeletonTwoSided(); }
80 bool requireUVVolume()
const
81 {
return local_assembler.doAlphaVolume(); }
82 bool requireUVSkeleton()
const
83 {
return local_assembler.doAlphaSkeleton(); }
84 bool requireUVBoundary()
const
85 {
return local_assembler.doAlphaBoundary(); }
86 bool requireUVVolumePostSkeleton()
const
87 {
return local_assembler.doAlphaVolumePostSkeleton(); }
93 return local_assembler;
97 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints&
trialConstraints()
const
103 const typename LocalAssembler::Traits::TestGridFunctionSpaceConstraints&
testConstraints()
const
112 global_a_ss_view.attach(jacobian_);
113 global_a_sn_view.attach(jacobian_);
114 global_a_ns_view.attach(jacobian_);
115 global_a_nn_view.attach(jacobian_);
122 global_s_s_view.attach(solution_);
123 global_s_n_view.attach(solution_);
129 template<
typename EG,
typename LFSUC,
typename LFSVC>
130 void onBindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
132 global_s_s_view.bind(lfsu_cache);
133 xl.
resize(lfsu_cache.size());
134 global_a_ss_view.bind(lfsv_cache,lfsu_cache);
135 al.assign(lfsv_cache.size(),lfsu_cache.size(),0.0);
138 template<
typename IG,
typename LFSUC,
typename LFSVC>
139 void onBindLFSUVOutside(
const IG & ig,
140 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
141 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
143 global_s_n_view.bind(lfsu_n_cache);
144 xn.
resize(lfsu_n_cache.size());
145 global_a_sn_view.bind(lfsv_s_cache,lfsu_n_cache);
146 al_sn.assign(lfsv_s_cache.size(),lfsu_n_cache.size(),0.0);
147 global_a_ns_view.bind(lfsv_n_cache,lfsu_s_cache);
148 al_ns.assign(lfsv_n_cache.size(),lfsu_s_cache.size(),0.0);
149 global_a_nn_view.bind(lfsv_n_cache,lfsu_n_cache);
150 al_nn.assign(lfsv_n_cache.size(),lfsu_n_cache.size(),0.0);
158 template<
typename EG,
typename LFSUC,
typename LFSVC>
159 void onUnbindLFSUV(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
161 local_assembler.scatter_jacobian(al,global_a_ss_view,
false);
164 template<
typename IG,
typename LFSUC,
typename LFSVC>
165 void onUnbindLFSUVOutside(
const IG & ig,
166 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
167 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
169 local_assembler.scatter_jacobian(al_sn,global_a_sn_view,
false);
170 local_assembler.scatter_jacobian(al_ns,global_a_ns_view,
false);
171 local_assembler.scatter_jacobian(al_nn,global_a_nn_view,
false);
178 template<
typename LFSUC>
181 global_s_s_view.read(xl);
183 template<
typename LFSUC>
184 void loadCoefficientsLFSUOutside(
const LFSUC & lfsu_n_cache)
186 global_s_n_view.read(xn);
188 template<
typename LFSUC>
189 void loadCoefficientsLFSUCoupling(
const LFSUC & lfsu_c_cache)
199 Jacobian& jacobian = global_a_ss_view.container();
200 global_s_s_view.detach();
201 global_s_n_view.detach();
202 global_a_ss_view.detach();
203 global_a_sn_view.detach();
204 global_a_ns_view.detach();
205 global_a_nn_view.detach();
207 if(local_assembler.doPostProcessing())
208 local_assembler.handle_dirichlet_constraints(gfsv,jacobian);
222 template<
typename EG>
226 (LocalAssembler::isNonOverlapping &&
236 template<
typename IG>
242 template<
typename EG,
typename LFSUC,
typename LFSVC>
243 void assembleUVVolume(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
245 al_view.setWeight(local_assembler.weight());
246 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolume>::
247 jacobian_volume(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),al_view);
250 template<
typename IG,
typename LFSUC,
typename LFSVC>
251 void assembleUVSkeleton(
const IG & ig,
const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
252 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache)
254 al_view.setWeight(local_assembler.weight());
255 al_sn_view.setWeight(local_assembler.weight());
256 al_ns_view.setWeight(local_assembler.weight());
257 al_nn_view.setWeight(local_assembler.weight());
259 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaSkeleton>::
260 jacobian_skeleton(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),lfsu_n_cache.localFunctionSpace(),xn,lfsv_n_cache.localFunctionSpace(),al_view,al_sn_view,al_ns_view,al_nn_view);
263 template<
typename IG,
typename LFSUC,
typename LFSVC>
264 void assembleUVBoundary(
const IG & ig,
const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache)
266 al_view.setWeight(local_assembler.weight());
267 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaBoundary>::
268 jacobian_boundary(lop,ig,lfsu_s_cache.localFunctionSpace(),xl,lfsv_s_cache.localFunctionSpace(),al_view);
271 template<
typename IG,
typename LFSUC,
typename LFSVC>
272 static void assembleUVEnrichedCoupling(
const IG & ig,
273 const LFSUC & lfsu_s_cache,
const LFSVC & lfsv_s_cache,
274 const LFSUC & lfsu_n_cache,
const LFSVC & lfsv_n_cache,
275 const LFSUC & lfsu_coupling_cache,
const LFSVC & lfsv_coupling_cache)
280 template<
typename IG,
typename LFSVC>
281 static void assembleVEnrichedCoupling(
const IG & ig,
282 const LFSVC & lfsv_s_cache,
283 const LFSVC & lfsv_n_cache,
284 const LFSVC & lfsv_coupling_cache)
289 template<
typename EG,
typename LFSUC,
typename LFSVC>
290 void assembleUVVolumePostSkeleton(
const EG & eg,
const LFSUC & lfsu_cache,
const LFSVC & lfsv_cache)
292 al_view.setWeight(local_assembler.weight());
293 Dune::PDELab::LocalAssemblerCallSwitch<LOP,LOP::doAlphaVolumePostSkeleton>::
294 jacobian_volume_post_skeleton(lop,eg,lfsu_cache.localFunctionSpace(),xl,lfsv_cache.localFunctionSpace(),al_view);
308 SolutionView global_s_s_view;
309 SolutionView global_s_n_view;
312 JacobianView global_a_ss_view;
313 JacobianView global_a_sn_view;
314 JacobianView global_a_ns_view;
315 JacobianView global_a_nn_view;
323 typedef typename std::conditional<
325 lop::DiagonalJacobian,
330 >::type JacobianMatrix;
336 JacobianMatrix al_sn;
337 JacobianMatrix al_ns;
338 JacobianMatrix al_nn;
340 typename JacobianMatrix::WeightedAccumulationView al_view;
341 typename JacobianMatrix::WeightedAccumulationView al_sn_view;
342 typename JacobianMatrix::WeightedAccumulationView al_ns_view;
343 typename JacobianMatrix::WeightedAccumulationView al_nn_view;
Default exception for dummy implementations.
Definition: exceptions.hh:263
The local assembler engine for DUNE grids which assembles the jacobian matrix.
Definition: jacobianengine.hh:26
void loadCoefficientsLFSUInside(const LFSUC &lfsu_cache)
Definition: jacobianengine.hh:179
void onUnbindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianengine.hh:159
const LocalAssembler & localAssembler() const
Public access to the wrapping local assembler.
Definition: jacobianengine.hh:91
bool requireSkeleton() const
Definition: jacobianengine.hh:76
void setSolution(const Solution &solution_)
Definition: jacobianengine.hh:120
void postAssembly(const GFSU &gfsu, const GFSV &gfsv)
Definition: jacobianengine.hh:197
bool skipIntersection(const IG &ig)
Definition: jacobianengine.hh:237
LA LocalAssembler
The type of the wrapping local assembler.
Definition: jacobianengine.hh:36
void onBindLFSUV(const EG &eg, const LFSUC &lfsu_cache, const LFSVC &lfsv_cache)
Definition: jacobianengine.hh:130
LA::Traits::Solution Solution
The type of the solution vector.
Definition: jacobianengine.hh:55
void setJacobian(Jacobian &jacobian_)
Definition: jacobianengine.hh:110
bool skipEntity(const EG &eg)
Definition: jacobianengine.hh:223
LA::LocalOperator LOP
The type of the local operator.
Definition: jacobianengine.hh:39
const LocalAssembler::Traits::TestGridFunctionSpaceConstraints & testConstraints() const
Test space constraints.
Definition: jacobianengine.hh:103
LA::LFSU LFSU
The local function spaces.
Definition: jacobianengine.hh:42
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints & trialConstraints() const
Trial space constraints.
Definition: jacobianengine.hh:97
DefaultLocalJacobianAssemblerEngine(const LocalAssembler &local_assembler_)
Constructor.
Definition: jacobianengine.hh:65
LA::Traits::Jacobian Jacobian
The type of the jacobian matrix.
Definition: jacobianengine.hh:50
A dense matrix for storing data associated with the degrees of freedom of a pair of LocalFunctionSpac...
Definition: diagonallocalmatrix.hh:29
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
void resize(size_type size)
Resize the container.
Definition: localvector.hh:324
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
Dune namespace.
Definition: alignedallocator.hh:13
Definition: localfunctionspacetags.hh:54
Definition: localfunctionspacetags.hh:48