1#ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
2#define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
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>
21 template<
typename GFSU,
typename GFSV,
typename CU,
typename CV>
28 using Element =
typename EntitySet::Element;
29 using Intersection =
typename EntitySet::Intersection;
35 typedef GFSV TestGridFunctionSpace;
39 typedef typename GFSU::Traits::SizeType
SizeType;
44 DefaultAssembler (
const GFSU& gfsu_,
const GFSV& gfsv_,
const CU& cu_,
const CV& cv_)
84 template<
class LocalAssemblerEngine>
87 typedef LFSIndexCache<LFSU,CU> LFSUCache;
89 typedef LFSIndexCache<LFSV,CV> LFSVCache;
91 const bool needs_constraints_caching = assembler_engine.needsConstraintsCaching(cu,cv);
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);
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();
112 auto entity_set = gfsu.entitySet();
113 auto& index_set = entity_set.indexSet();
116 for (
const auto& element : elements(entity_set))
119 auto ids = index_set.uniqueIndex(element);
127 lfsv.bind( element );
131 assembler_engine.onBindLFSV(eg,lfsv_cache);
137 lfsu.bind( element );
141 assembler_engine.onBindLFSUV(eg,lfsu_cache,lfsv_cache);
144 assembler_engine.loadCoefficientsLFSUInside(lfsu_cache);
150 if (require_uv_skeleton || require_v_skeleton ||
151 require_uv_boundary || require_v_boundary ||
152 require_uv_processor || require_v_processor)
155 unsigned int intersection_index = 0;
156 for(
const auto& intersection : intersections(entity_set,element))
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);
168 switch (intersection_type)
170 case IntersectionType::skeleton:
173 case IntersectionType::periodic:
174 if (require_uv_skeleton || require_v_skeleton)
178 auto idn = index_set.uniqueIndex(outside_element);
181 bool visit_face = ids > idn || require_skeleton_two_sided;
187 lfsvn.bind(outside_element);
188 lfsvn_cache.update();
191 assembler_engine.onBindLFSVOutside(ig,lfsv_cache,lfsvn_cache);
196 if(require_uv_skeleton){
199 lfsun.bind(outside_element);
200 lfsun_cache.update();
203 assembler_engine.onBindLFSUVOutside(ig,
204 lfsu_cache,lfsv_cache,
205 lfsun_cache,lfsvn_cache);
208 assembler_engine.loadCoefficientsLFSUOutside(lfsun_cache);
211 assembler_engine.
assembleUVSkeleton(ig,lfsu_cache,lfsv_cache,lfsun_cache,lfsvn_cache);
214 assembler_engine.onUnbindLFSUVOutside(ig,
215 lfsu_cache,lfsv_cache,
216 lfsun_cache,lfsvn_cache);
220 assembler_engine.onUnbindLFSVOutside(ig,lfsv_cache,lfsvn_cache);
225 case IntersectionType::boundary:
226 if(require_uv_boundary || require_v_boundary )
232 if(require_uv_boundary){
239 case IntersectionType::processor:
240 if(require_uv_processor || require_v_processor )
246 if(require_uv_processor){
254 ++intersection_index;
258 if(require_uv_post_skeleton || require_v_post_skeleton){
262 if(require_uv_post_skeleton){
269 assembler_engine.onUnbindLFSUV(eg,lfsu_cache,lfsv_cache);
272 assembler_engine.onUnbindLFSV(eg,lfsv_cache);
287 typename std::conditional<
288 std::is_same<CU,EmptyTransformation>::value,
292 typename std::conditional<
293 std::is_same<CV,EmptyTransformation>::value,
299 typedef LocalFunctionSpace<GFSU, TrialSpaceTag> LFSU;
300 typedef LocalFunctionSpace<GFSV, TestSpaceTag> LFSV;
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)
bool skipEntity(const EG &eg)
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:13