DUNE PDELab (2.8)

assembler.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
3
5#include <dune/pdelab/gridoperator/common/assemblerutilities.hh>
6#include <dune/pdelab/gridfunctionspace/localfunctionspace.hh>
7#include <dune/pdelab/gridfunctionspace/lfsindexcache.hh>
8#include <dune/pdelab/common/geometrywrapper.hh>
9#include <dune/pdelab/common/intersectiontype.hh>
10
11namespace Dune{
12 namespace PDELab{
13
21 template<typename GFSU, typename GFSV, typename CU, typename CV>
23 public:
24
27 using EntitySet = typename GFSU::Traits::EntitySet;
28 using Element = typename EntitySet::Element;
29 using Intersection = typename EntitySet::Intersection;
31
35 typedef GFSV TestGridFunctionSpace;
37
39 typedef typename GFSU::Traits::SizeType SizeType;
40
42 static const bool isGalerkinMethod = std::is_same<GFSU,GFSV>::value;
43
44 DefaultAssembler (const GFSU& gfsu_, const GFSV& gfsv_, const CU& cu_, const CV& cv_)
45 : gfsu(gfsu_)
46 , gfsv(gfsv_)
47 , cu(cu_)
48 , cv(cv_)
49 , lfsu(gfsu_)
50 , lfsv(gfsv_)
51 , lfsun(gfsu_)
52 , lfsvn(gfsv_)
53 { }
54
55 DefaultAssembler (const GFSU& gfsu_, const GFSV& gfsv_)
56 : gfsu(gfsu_)
57 , gfsv(gfsv_)
58 , cu()
59 , cv()
60 , lfsu(gfsu_)
61 , lfsv(gfsv_)
62 , lfsun(gfsu_)
63 , lfsvn(gfsv_)
64 { }
65
67 const GFSU& trialGridFunctionSpace() const
68 {
69 return gfsu;
70 }
71
73 const GFSV& testGridFunctionSpace() const
74 {
75 return gfsv;
76 }
77
78 // Assembler (const GFSU& gfsu_, const GFSV& gfsv_)
79 // : gfsu(gfsu_), gfsv(gfsv_), lfsu(gfsu_), lfsv(gfsv_),
80 // lfsun(gfsu_), lfsvn(gfsv_),
81 // sub_triangulation(ST(gfsu_.gridview(),Dune::PDELab::NoSubTriangulationImp()))
82 // { }
83
84 template<class LocalAssemblerEngine>
85 void assemble(LocalAssemblerEngine & assembler_engine) const
86 {
87 typedef LFSIndexCache<LFSU,CU> LFSUCache;
88
89 typedef LFSIndexCache<LFSV,CV> LFSVCache;
90
91 const bool needs_constraints_caching = assembler_engine.needsConstraintsCaching(cu,cv);
92
93 LFSUCache lfsu_cache(lfsu,cu,needs_constraints_caching);
94 LFSVCache lfsv_cache(lfsv,cv,needs_constraints_caching);
95 LFSUCache lfsun_cache(lfsun,cu,needs_constraints_caching);
96 LFSVCache lfsvn_cache(lfsvn,cv,needs_constraints_caching);
97
98 // Notify assembler engine about oncoming assembly
99 assembler_engine.preAssembly();
100
101 // Extract integration requirements from the local assembler
102 const bool require_uv_skeleton = assembler_engine.requireUVSkeleton();
103 const bool require_v_skeleton = assembler_engine.requireVSkeleton();
104 const bool require_uv_boundary = assembler_engine.requireUVBoundary();
105 const bool require_v_boundary = assembler_engine.requireVBoundary();
106 const bool require_uv_processor = assembler_engine.requireUVBoundary();
107 const bool require_v_processor = assembler_engine.requireVBoundary();
108 const bool require_uv_post_skeleton = assembler_engine.requireUVVolumePostSkeleton();
109 const bool require_v_post_skeleton = assembler_engine.requireVVolumePostSkeleton();
110 const bool require_skeleton_two_sided = assembler_engine.requireSkeletonTwoSided();
111
112 auto entity_set = gfsu.entitySet();
113 auto& index_set = entity_set.indexSet();
114
115 // Traverse grid view
116 for (const auto& element : elements(entity_set))
117 {
118 // Compute unique id
119 auto ids = index_set.uniqueIndex(element);
120
121 ElementGeometry<Element> eg(element);
122
123 if(assembler_engine.skipEntity(eg))
124 continue;
125
126 // Bind local test function space to element
127 lfsv.bind( element );
128 lfsv_cache.update();
129
130 // Notify assembler engine about bind
131 assembler_engine.onBindLFSV(eg,lfsv_cache);
132
133 // Volume integration
134 assembler_engine.assembleVVolume(eg,lfsv_cache);
135
136 // Bind local trial function space to element
137 lfsu.bind( element );
138 lfsu_cache.update();
139
140 // Notify assembler engine about bind
141 assembler_engine.onBindLFSUV(eg,lfsu_cache,lfsv_cache);
142
143 // Load coefficients of local functions
144 assembler_engine.loadCoefficientsLFSUInside(lfsu_cache);
145
146 // Volume integration
147 assembler_engine.assembleUVVolume(eg,lfsu_cache,lfsv_cache);
148
149 // Skip if no intersection iterator is needed
150 if (require_uv_skeleton || require_v_skeleton ||
151 require_uv_boundary || require_v_boundary ||
152 require_uv_processor || require_v_processor)
153 {
154 // Traverse intersections
155 unsigned int intersection_index = 0;
156 for(const auto& intersection : intersections(entity_set,element))
157 {
158
159 IntersectionGeometry<Intersection> ig(intersection,intersection_index);
160
161 if(assembler_engine.skipIntersection(ig))
162 continue;
163
164 auto intersection_data = classifyIntersection(entity_set,intersection);
165 auto intersection_type = std::get<0>(intersection_data);
166 auto& outside_element = std::get<1>(intersection_data);
167
168 switch (intersection_type)
169 {
170 case IntersectionType::skeleton:
171 // the specific ordering of the if-statements in the old code caused periodic
172 // boundary intersection to be handled the same as skeleton intersections
173 case IntersectionType::periodic:
174 if (require_uv_skeleton || require_v_skeleton)
175 {
176 // compute unique id for neighbor
177
178 auto idn = index_set.uniqueIndex(outside_element);
179
180 // Visit face if id is bigger
181 bool visit_face = ids > idn || require_skeleton_two_sided;
182
183 // unique vist of intersection
184 if (visit_face)
185 {
186 // Bind local test space to neighbor element
187 lfsvn.bind(outside_element);
188 lfsvn_cache.update();
189
190 // Notify assembler engine about binds
191 assembler_engine.onBindLFSVOutside(ig,lfsv_cache,lfsvn_cache);
192
193 // Skeleton integration
194 assembler_engine.assembleVSkeleton(ig,lfsv_cache,lfsvn_cache);
195
196 if(require_uv_skeleton){
197
198 // Bind local trial space to neighbor element
199 lfsun.bind(outside_element);
200 lfsun_cache.update();
201
202 // Notify assembler engine about binds
203 assembler_engine.onBindLFSUVOutside(ig,
204 lfsu_cache,lfsv_cache,
205 lfsun_cache,lfsvn_cache);
206
207 // Load coefficients of local functions
208 assembler_engine.loadCoefficientsLFSUOutside(lfsun_cache);
209
210 // Skeleton integration
211 assembler_engine.assembleUVSkeleton(ig,lfsu_cache,lfsv_cache,lfsun_cache,lfsvn_cache);
212
213 // Notify assembler engine about unbinds
214 assembler_engine.onUnbindLFSUVOutside(ig,
215 lfsu_cache,lfsv_cache,
216 lfsun_cache,lfsvn_cache);
217 }
218
219 // Notify assembler engine about unbinds
220 assembler_engine.onUnbindLFSVOutside(ig,lfsv_cache,lfsvn_cache);
221 }
222 }
223 break;
224
225 case IntersectionType::boundary:
226 if(require_uv_boundary || require_v_boundary )
227 {
228
229 // Boundary integration
230 assembler_engine.assembleVBoundary(ig,lfsv_cache);
231
232 if(require_uv_boundary){
233 // Boundary integration
234 assembler_engine.assembleUVBoundary(ig,lfsu_cache,lfsv_cache);
235 }
236 }
237 break;
238
239 case IntersectionType::processor:
240 if(require_uv_processor || require_v_processor )
241 {
242
243 // Processor integration
244 assembler_engine.assembleVProcessor(ig,lfsv_cache);
245
246 if(require_uv_processor){
247 // Processor integration
248 assembler_engine.assembleUVProcessor(ig,lfsu_cache,lfsv_cache);
249 }
250 }
251 break;
252 } // switch
253
254 ++intersection_index;
255 } // iit
256 } // do skeleton
257
258 if(require_uv_post_skeleton || require_v_post_skeleton){
259 // Volume integration
260 assembler_engine.assembleVVolumePostSkeleton(eg,lfsv_cache);
261
262 if(require_uv_post_skeleton){
263 // Volume integration
264 assembler_engine.assembleUVVolumePostSkeleton(eg,lfsu_cache,lfsv_cache);
265 }
266 }
267
268 // Notify assembler engine about unbinds
269 assembler_engine.onUnbindLFSUV(eg,lfsu_cache,lfsv_cache);
270
271 // Notify assembler engine about unbinds
272 assembler_engine.onUnbindLFSV(eg,lfsv_cache);
273
274 } // it
275
276 // Notify assembler engine that assembly is finished
277 assembler_engine.postAssembly(gfsu,gfsv);
278
279 }
280
281 private:
282
283 /* global function spaces */
284 const GFSU& gfsu;
285 const GFSV& gfsv;
286
287 typename std::conditional<
288 std::is_same<CU,EmptyTransformation>::value,
289 const CU,
290 const CU&
291 >::type cu;
292 typename std::conditional<
293 std::is_same<CV,EmptyTransformation>::value,
294 const CV,
295 const CV&
296 >::type cv;
297
298 /* local function spaces */
299 typedef LocalFunctionSpace<GFSU, TrialSpaceTag> LFSU;
300 typedef LocalFunctionSpace<GFSV, TestSpaceTag> LFSV;
301 // local function spaces in local cell
302 mutable LFSU lfsu;
303 mutable LFSV lfsv;
304 // local function spaces in neighbor
305 mutable LFSU lfsun;
306 mutable LFSV lfsvn;
307
308 };
309
310 }
311}
312#endif // DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
The assembler for standard DUNE grid.
Definition: assembler.hh:22
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: assembler.hh:67
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: assembler.hh:73
GFSU::Traits::SizeType SizeType
Size type as used in grid function space.
Definition: assembler.hh:39
typename GFSU::Traits::EntitySet EntitySet
Definition: assembler.hh:27
GFSU TrialGridFunctionSpace
Definition: assembler.hh:34
static const bool isGalerkinMethod
Static check on whether this is a Galerkin method.
Definition: assembler.hh:42
Wrap element.
Definition: geometrywrapper.hh:16
Wrap intersection.
Definition: geometrywrapper.hh:57
The local assembler engine which handles the integration parts as provided by the global assemblers.
Definition: assembler.hh:34
void assembleUVVolume(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
void assembleVSkeleton(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void assembleUVProcessor(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
bool skipIntersection(const IG &ig)
void assembleUVSkeleton(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n)
void assembleVProcessor(const IG &ig, const LFSV_S &lfsv_s)
void preAssembly()
Called directly before assembling.
void assembleVVolumePostSkeleton(const EG &eg, const LFSV &lfsv)
void assembleUVVolumePostSkeleton(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
void assembleVVolume(const EG &eg, const LFSV &lfsv)
void postAssembly()
Called last thing after assembling.
void assembleUVBoundary(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
Traits for type conversions and type information.
Dune namespace.
Definition: alignedallocator.hh:11
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)