1 #ifndef DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_DEFAULT_ASSEMBLER_HH
4 #include <dune/common/typetraits.hh>
21 template<
typename GFSU,
typename GFSV,
typename CU,
typename CV>
28 using Element =
typename EntitySet::Element;
39 typedef typename GFSU::Traits::SizeType
SizeType;
84 template<
class LocalAssemblerEngine>
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);
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 );
137 lfsu.bind( element );
141 assembler_engine.
onBindLFSUV(eg,lfsu_cache,lfsv_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))
165 auto intersection_type = std::get<0>(intersection_data);
166 auto& outside_element = std::get<1>(intersection_data);
168 switch (intersection_type)
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();
196 if(require_uv_skeleton){
199 lfsun.bind(outside_element);
200 lfsun_cache.update();
204 lfsu_cache,lfsv_cache,
205 lfsun_cache,lfsvn_cache);
215 lfsu_cache,lfsv_cache,
216 lfsun_cache,lfsvn_cache);
226 if(require_uv_boundary || require_v_boundary )
232 if(require_uv_boundary){
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){
287 typename std::conditional<
292 typename std::conditional<
const IG & ig
Definition: constraints.hh:149
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
@ processor
processor boundary intersection (neighbor() == false && boundary() == false) or outside entity not in...
@ boundary
domain boundary intersection (neighbor() == false && boundary() == true)
@ skeleton
skeleton intersection (neighbor() == true && boundary() == false)
@ periodic
periodic boundary intersection (neighbor() == true && boundary() == true)
std::tuple< IntersectionType, typename EntitySet::Element > classifyIntersection(const EntitySet &entity_set, const Intersection &is)
Classifies the type of an intersection wrt to the passed EntitySet.
Definition: intersectiontype.hh:37
Wrap element.
Definition: geometrywrapper.hh:16
Wrap intersection.
Definition: geometrywrapper.hh:57
Definition: lfsindexcache.hh:979
The local assembler engine which handles the integration parts as provided by the global assemblers.
Definition: common/assembler.hh:34
void onUnbindLFSV(const EG &eg, const LFSV_S &lfsv_s)
void assembleUVVolume(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
bool requireSkeletonTwoSided() const
void onBindLFSUV(const EG &eg, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void assembleVSkeleton(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
bool requireVSkeleton() const
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)
bool requireUVVolumePostSkeleton() const
bool requireVBoundary() const
void loadCoefficientsLFSUInside(const LFSU_S &lfsu_s)
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 onUnbindLFSUVOutside(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n)
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
void onBindLFSV(const EG &eg, const LFSV_S &lfsv_s)
void onBindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void onUnbindLFSUV(const EG &eg, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
bool requireUVSkeleton() const
void assembleVVolume(const EG &eg, const LFSV &lfsv)
void onUnbindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void postAssembly()
Called last thing after assembling.
bool requireUVBoundary() const
void assembleUVBoundary(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void onBindLFSUVOutside(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n)
void loadCoefficientsLFSUOutside(const LFSU_N &lfsu_n)
bool requireVVolumePostSkeleton() const
The assembler for standard DUNE grid.
Definition: default/assembler.hh:22
void assemble(LocalAssemblerEngine &assembler_engine) const
Definition: default/assembler.hh:85
typename EntitySet::Element Element
Definition: default/assembler.hh:28
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: default/assembler.hh:67
GFSU::Traits::SizeType SizeType
Size type as used in grid function space.
Definition: default/assembler.hh:39
DefaultAssembler(const GFSU &gfsu_, const GFSV &gfsv_, const CU &cu_, const CV &cv_)
Definition: default/assembler.hh:44
GFSV TestGridFunctionSpace
Definition: default/assembler.hh:35
typename GFSU::Traits::EntitySet EntitySet
Definition: default/assembler.hh:27
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: default/assembler.hh:73
DefaultAssembler(const GFSU &gfsu_, const GFSV &gfsv_)
Definition: default/assembler.hh:55
typename EntitySet::Intersection Intersection
Definition: default/assembler.hh:29
GFSU TrialGridFunctionSpace
Definition: default/assembler.hh:34
static const bool isGalerkinMethod
Static check on whether this is a Galerkin method.
Definition: default/assembler.hh:42
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139