Go to the documentation of this file.
3 #ifndef DUNE_PDELAB_GRIDOPERATOR_FASTDG_ASSEMBLER_HH
4 #define DUNE_PDELAB_GRIDOPERATOR_FASTDG_ASSEMBLER_HH
6 #include <dune/common/typetraits.hh>
24 template<
typename GFSU,
typename GFSV,
typename CU,
typename CV,
bool nonoverlapping_mode=false>
31 using Element =
typename EntitySet::Element;
42 typedef typename GFSU::Traits::SizeType
SizeType;
47 FastDGAssembler (
const GFSU& gfsu_,
const GFSV& gfsv_,
const CU& cu_,
const CV& cv_)
92 template<
class EngineFactory,
class LocalAssembler>
93 void assemble(
const EngineFactory &engineFactory,
94 LocalAssembler &la)
const
99 template<
class LocalAssemblerEngine>
102 const bool fast =
true;
107 const bool needs_constraints_caching = assembler_engine.needsConstraintsCaching(cu,cv);
109 LFSUCache lfsu_cache(lfsu,cu,needs_constraints_caching);
110 LFSVCache lfsv_cache(lfsv,cv,needs_constraints_caching);
111 LFSUCache lfsun_cache(lfsun,cu,needs_constraints_caching);
112 LFSVCache lfsvn_cache(lfsvn,cv,needs_constraints_caching);
128 auto entity_set = gfsu.entitySet();
129 auto& index_set = entity_set.indexSet();
132 for (
const auto& element : elements(entity_set,assembler_engine.partition()))
135 auto ids = index_set.uniqueIndex(element);
143 lfsv.bind(element,std::integral_constant<bool,fast>{});
153 lfsu.bind(element,std::integral_constant<bool,fast>{});
157 assembler_engine.
onBindLFSUV(eg,lfsu_cache,lfsv_cache);
166 if (require_uv_skeleton || require_v_skeleton ||
167 require_uv_boundary || require_v_boundary ||
168 require_uv_processor || require_v_processor)
171 unsigned int intersection_index = 0;
172 for(
const auto& intersection : intersections(entity_set,element))
178 auto intersection_type = std::get<0>(intersection_data);
179 auto& outside_element = std::get<1>(intersection_data);
181 switch (intersection_type)
187 if (require_uv_skeleton || require_v_skeleton)
191 auto idn = index_set.uniqueIndex(outside_element);
194 bool visit_face = ids > idn
195 or require_skeleton_two_sided
196 or (assembler_engine.partition() == Partitions::interiorBorder
197 and outside_element.partitionType() != InteriorEntity);
203 lfsvn.bind(outside_element,std::integral_constant<bool,fast>{});
204 lfsvn_cache.update();
212 if(require_uv_skeleton){
215 lfsun.bind(outside_element,std::integral_constant<bool,fast>{});
216 lfsun_cache.update();
220 lfsu_cache,lfsv_cache,
221 lfsun_cache,lfsvn_cache);
231 lfsu_cache,lfsv_cache,
232 lfsun_cache,lfsvn_cache);
242 if(require_uv_boundary || require_v_boundary )
248 if(require_uv_boundary){
256 if(require_uv_processor || require_v_processor )
262 if(require_uv_processor){
270 ++intersection_index;
274 if(require_uv_post_skeleton || require_v_post_skeleton){
278 if(require_uv_post_skeleton){
303 typename std::conditional<
308 typename std::conditional<
bool requireSkeletonTwoSided() const
void onBindLFSUV(const EG &eg, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void onUnbindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
typename GFSU::Traits::EntitySet EntitySet
Definition: fastdg/assembler.hh:30
GFSU TrialGridFunctionSpace
Definition: fastdg/assembler.hh:37
void preAssembly()
Called directly before assembling.
typename EntitySet::Intersection Intersection
Definition: fastdg/assembler.hh:32
typename EntitySet::Element Element
Definition: fastdg/assembler.hh:31
GFSU::Traits::SizeType SizeType
Size type as used in grid function space.
Definition: fastdg/assembler.hh:42
void assembleUVVolumePostSkeleton(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
bool requireVSkeleton() const
FastDGAssembler(const GFSU &gfsu_, const GFSV &gfsv_, const CU &cu_, const CV &cv_)
Definition: fastdg/assembler.hh:47
bool requireUVSkeleton() const
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
void loadCoefficientsLFSUOutside(const LFSU_N &lfsu_n)
Wrap element.
Definition: geometrywrapper.hh:15
void postAssembly()
Called last thing after assembling.
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 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 assembleUVProcessor(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void onUnbindLFSUV(const EG &eg, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
static const bool isGalerkinMethod
Static check on whether this is a Galerkin method.
Definition: fastdg/assembler.hh:45
The fast DG assembler for standard DUNE grid.
Definition: fastdg/assembler.hh:25
void loadCoefficientsLFSUInside(const LFSU_S &lfsu_s)
bool requireUVVolumePostSkeleton() const
const GFSV & testGridFunctionSpace() const
Get the test grid function space.
Definition: fastdg/assembler.hh:76
void assembleVVolume(const EG &eg, const LFSV &lfsv)
const IG & ig
Definition: constraints.hh:148
void assembleVVolumePostSkeleton(const EG &eg, const LFSV &lfsv)
@ processor
processor boundary intersection (neighbor() == false && boundary() == false) or outside entity not in...
void assemble(LocalAssemblerEngine &assembler_engine) const
Definition: fastdg/assembler.hh:100
@ skeleton
skeleton intersection (neighbor() == true && boundary() == false)
The local assembler engine which handles the integration parts as provided by the global assemblers.
Definition: common/assembler.hh:33
@ periodic
periodic boundary intersection (neighbor() == true && boundary() == true)
Wrap intersection.
Definition: geometrywrapper.hh:56
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 assembleUVBoundary(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
void assemble(const EngineFactory &engineFactory, LocalAssembler &la) const
do the assembly
Definition: fastdg/assembler.hh:93
@ boundary
domain boundary intersection (neighbor() == false && boundary() == true)
void assembleVProcessor(const IG &ig, const LFSV_S &lfsv_s)
void onUnbindLFSV(const EG &eg, const LFSV_S &lfsv_s)
bool requireUVBoundary() const
bool requireVBoundary() const
const GFSU & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: fastdg/assembler.hh:70
GFSV TestGridFunctionSpace
Definition: fastdg/assembler.hh:38
bool requireVVolumePostSkeleton() const
FastDGAssembler(const GFSU &gfsu_, const GFSV &gfsv_)
Definition: fastdg/assembler.hh:58
bool assembleCell(const EG &eg)
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
void assembleUVVolume(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
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 assembleVSkeleton(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
void onBindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void onBindLFSV(const EG &eg, const LFSV_S &lfsv_s)
Definition: lfsindexcache.hh:977