1 #ifndef DUNE_PDELAB_DEFAULT_ASSEMBLER_HH 2 #define DUNE_PDELAB_DEFAULT_ASSEMBLER_HH 4 #include <dune/common/typetraits.hh> 22 template<
typename GFSU,
typename GFSV,
typename CU,
typename CV,
bool nonoverlapping_mode=false>
29 using Element =
typename EntitySet::Element;
40 typedef typename GFSU::Traits::SizeType
SizeType;
85 template<
class LocalAssemblerEngine>
92 const bool needs_constraints_caching = assembler_engine.needsConstraintsCaching(cu,cv);
94 LFSUCache lfsu_cache(lfsu,cu,needs_constraints_caching);
95 LFSVCache lfsv_cache(lfsv,cv,needs_constraints_caching);
96 LFSUCache lfsun_cache(lfsun,cu,needs_constraints_caching);
97 LFSVCache lfsvn_cache(lfsvn,cv,needs_constraints_caching);
113 auto entity_set = gfsu.entitySet();
114 auto& index_set = entity_set.indexSet();
117 for (
const auto& element : elements(entity_set))
120 auto ids = index_set.uniqueIndex(element);
128 lfsv.bind( element );
138 lfsu.bind( element );
142 assembler_engine.
onBindLFSUV(eg,lfsu_cache,lfsv_cache);
151 if (require_uv_skeleton || require_v_skeleton ||
152 require_uv_boundary || require_v_boundary ||
153 require_uv_processor || require_v_processor)
156 unsigned int intersection_index = 0;
157 for(
const auto& intersection : intersections(entity_set,element))
163 auto intersection_type = std::get<0>(intersection_data);
164 auto& outside_element = std::get<1>(intersection_data);
166 switch (intersection_type)
172 if (require_uv_skeleton || require_v_skeleton)
176 auto idn = index_set.uniqueIndex(outside_element);
179 bool visit_face = ids > idn || require_skeleton_two_sided;
185 lfsvn.bind(outside_element);
186 lfsvn_cache.update();
194 if(require_uv_skeleton){
197 lfsun.bind(outside_element);
198 lfsun_cache.update();
202 lfsu_cache,lfsv_cache,
203 lfsun_cache,lfsvn_cache);
209 assembler_engine.
assembleUVSkeleton(ig,lfsu_cache,lfsv_cache,lfsun_cache,lfsvn_cache);
213 lfsu_cache,lfsv_cache,
214 lfsun_cache,lfsvn_cache);
224 if(require_uv_boundary || require_v_boundary )
230 if(require_uv_boundary){
238 if(require_uv_processor || require_v_processor )
244 if(require_uv_processor){
252 ++intersection_index;
256 if(require_uv_post_skeleton || require_v_post_skeleton){
260 if(require_uv_post_skeleton){
285 typename conditional<
290 typename conditional<
void assembleUVProcessor(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
skeleton intersection (neighbor() == true && boundary() == false)
void onUnbindLFSV(const EG &eg, const LFSV_S &lfsv_s)
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: default/assembler.hh:74
bool requireUVSkeleton() const
static const unsigned int value
Definition: gridfunctionspace/tags.hh:177
GFSV TestGridFunctionSpace
Definition: default/assembler.hh:36
DefaultAssembler(const GFSU &gfsu_, const GFSV &gfsv_)
Definition: default/assembler.hh:56
domain boundary intersection (neighbor() == false && boundary() == true)
void assembleVVolume(const EG &eg, const LFSV &lfsv)
const IG & ig
Definition: constraints.hh:148
bool requireUVVolumePostSkeleton() const
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: default/assembler.hh:68
typename EntitySet::Element Element
Definition: default/assembler.hh:29
bool requireVSkeleton() const
void assembleVProcessor(const IG &ig, const LFSV_S &lfsv_s)
void onBindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
Definition: adaptivity.hh:27
GFSU TrialGridFunctionSpace
Definition: default/assembler.hh:35
void assembleUVVolumePostSkeleton(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
void assembleUVBoundary(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void onBindLFSUV(const EG &eg, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
bool assembleCell(const EG &eg)
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
periodic boundary intersection (neighbor() == true && boundary() == true)
typename EntitySet::Intersection Intersection
Definition: default/assembler.hh:30
The assembler for standard DUNE grid.
Definition: default/assembler.hh:23
void onUnbindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
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)
typename GFSU::Traits::EntitySet EntitySet
Definition: default/assembler.hh:28
void preAssembly()
Called directly before assembling.
bool requireSkeletonTwoSided() const
void onUnbindLFSUV(const EG &eg, 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)
Wrap element.
Definition: geometrywrapper.hh:15
Definition: lfsindexcache.hh:948
bool requireUVBoundary() const
bool requireVVolumePostSkeleton() const
Wrap intersection.
Definition: geometrywrapper.hh:56
void postAssembly()
Called last thing after assembling.
void loadCoefficientsLFSUInside(const LFSU_S &lfsu_s)
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
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 onBindLFSV(const EG &eg, const LFSV_S &lfsv_s)
bool requireVBoundary() const
GFSU::Traits::SizeType SizeType
Size type as used in grid function space.
Definition: default/assembler.hh:40
void assembleVVolumePostSkeleton(const EG &eg, const LFSV &lfsv)
void loadCoefficientsLFSUOutside(const LFSU_N &lfsu_n)
static const bool isGalerkinMethod
Static check on whether this is a Galerkin method.
Definition: default/assembler.hh:43
processor boundary intersection (neighbor() == false && boundary() == false) or outside entity not in...
DefaultAssembler(const GFSU &gfsu_, const GFSV &gfsv_, const CU &cu_, const CV &cv_)
Definition: default/assembler.hh:45
void assemble(LocalAssemblerEngine &assembler_engine) const
Definition: default/assembler.hh:86
void assembleUVVolume(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
The local assembler engine which handles the integration parts as provided by the global assemblers...
Definition: common/assembler.hh:31
void assembleVSkeleton(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)