1 #ifndef DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH 2 #define DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH 4 #include <assemblerutilties.hh> 24 template<
class LocalAssemblerEngine>
39 const LocalAssembler & localAssembler();
48 bool requireSkeleton()
const;
49 bool requireSkeletonTwoSided()
const;
50 bool requireUVVolume()
const;
51 bool requireVVolume()
const;
52 bool requireUVSkeleton()
const;
53 bool requireVSkeleton()
const;
54 bool requireUVBoundary()
const;
55 bool requireVBoundary()
const;
56 bool requireUVProcessor()
const;
57 bool requireVProcessor()
const;
58 bool requireUVEnrichedCoupling()
const;
59 bool requireVEnrichedCoupling()
const;
60 bool requireUVVolumePostSkeleton()
const;
61 bool requireVVolumePostSkeleton()
const;
84 bool assembleCell(
const EG & eg);
89 template<
typename EG,
typename LFSU,
typename LFSV>
90 void assembleUVVolume(
const EG & eg,
const LFSU & lfsu,
const LFSV & lfsv);
95 template<
typename EG,
typename LFSV>
96 void assembleVVolume(
const EG & eg,
const LFSV & lfsv);
101 template<
typename IG,
typename LFSU_S,
typename LFSV_S,
typename LFSU_N,
typename LFSV_N>
102 void assembleUVSkeleton(
const IG &
ig,
const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
103 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
108 template<
typename IG,
typename LFSV_S,
typename LFSV_N>
109 void assembleVSkeleton(
const IG &
ig,
const LFSV_S & lfsv_s,
const LFSV_N & lfsv_n);
114 template<
typename IG,
typename LFSU_S,
typename LFSV_S>
115 void assembleUVBoundary(
const IG &
ig,
const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
120 template<
typename IG,
typename LFSV_S>
121 void assembleVBoundary(
const IG &
ig,
const LFSV_S & lfsv_s);
128 template<
typename IG,
typename LFSU_S,
typename LFSV_S>
129 void assembleUVProcessor(
const IG &
ig,
const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
136 template<
typename IG,
typename LFSV_S>
137 void assembleVProcessor(
const IG &
ig,
const LFSV_S & lfsv_s);
139 template<
typename IG,
typename LFSU_S,
typename LFSV_S,
typename LFSU_N,
typename LFSV_N,
140 typename LFSU_C,
typename LFSV_C>
141 void assembleUVEnrichedCoupling(
const IG &
ig,
142 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
143 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n,
144 const LFSU_C & lfsu_c,
const LFSV_C & lfsv_c);
146 template<
typename IG,
typename LFSV_S,
typename LFSV_N,
typename LFSV_C>
147 void assembleVEnrichedCoupling(
const IG &
ig,
148 const LFSV_S & lfsv_s,
149 const LFSV_N & lfsv_n,
150 const LFSV_C & lfsv_c);
156 template<
typename EG,
typename LFSU,
typename LFSV>
157 void assembleUVVolumePostSkeleton(
const EG & eg,
const LFSU & lfsu,
const LFSV & lfsv);
163 template<
typename EG,
typename LFSV>
164 void assembleVVolumePostSkeleton(
const EG & eg,
const LFSV & lfsv);
184 void onBindLFSUV(
const EG & eg,
185 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
186 void onBindLFSV(
const EG & eg,
187 const LFSV_S & lfsv_s);
188 void onBindLFSUVInside(
const IG &
ig,
189 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
190 void onBindLFSVInside(
const IG &
ig,
191 const LFSV_S & lfsv_s);
192 void onBindLFSUVOutside(
const IG &
ig,
193 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
194 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
195 void onBindLFSVOutside(
const IG &
ig,
196 const LFSV_S & lfsv_s,
197 const LFSV_N & lfsv_n);
198 void onBindLFSUVCoupling(
const IG &
ig,
199 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
200 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n
201 const LFSU_Coupling & lfsu_coupling,
const LFSV_Coupling & lfsv_coupling);
202 void onBindLFSVCoupling(
const IG &
ig,
203 const LFSV_S & lfsv_s,
204 const LFSV_N & lfsv_n,
205 const LFSV_Coupling & lfsv_coupling);
207 void onUnbindLFSUV(
const EG & eg,
208 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
209 void onUnbindLFSV(
const EG & eg,
210 const LFSV_S & lfsv_s);
211 void onUnbindLFSUVInside(
const IG &
ig,
212 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
213 void onUnbindLFSVInside(
const IG &
ig,
214 const LFSV_S & lfsv_s);
215 void onUnbindLFSUVOutside(
const IG &
ig,
216 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
217 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
218 void onUnbindLFSVOutside(
const IG &
ig,
219 const LFSV_S & lfsv_s,
220 const LFSV_N & lfsv_n);
221 void onUnbindLFSUVCoupling(
const IG &
ig,
222 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
223 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n,
224 const LFSU_Coupling & lfsu_coupling,
const LFSV_Coupling & lfsv_coupling);
225 void onUnbindLFSVCoupling(
const IG &
ig,
226 const LFSV_S & lfsv_s,
227 const LFSV_N & lfsv_n,
228 const LFSV_Coupling & lfsv_coupling);
240 void loadCoefficientsLFSUInside(
const LFSU_S & lfsu_s);
241 void loadCoefficientsLFSUOutside(
const LFSU_N & lfsu_n);
242 void loadCoefficientsLFSUCoupling(
const LFSU_Coupling & lfsu_coupling);
255 void setSolution(
const X& x);
256 void setPattern(
const P&
p);
257 void setJacobian(
const J & j);
258 void setResidual(
const R& r);
274 template<
typename B,
typename CU,
typename CV>
285 void setTime(TT time);
288 template<
typename TT>
289 void preStep (TT time, TT dt, std::size_t stages);
295 template<
typename TT>
296 void preStage (TT time, std::size_t stage);
302 template<
typename TT>
303 TT suggestTimestep (TT dt)
const;
309 void setWeight(RF weight);
344 template<
typename GFSU,
typename GFSV,
345 typename MB,
typename DF,
typename RF,
typename JF>
355 void fill_pattern (P& globalpattern)
const;
359 template<
typename X,
typename R>
360 void residual (
const X& x, R& r)
const;
364 template<
typename X,
typename A>
365 void jacobian (
const X& x, A& a)
const;
369 Assembler & assembler();
375 const GFSU& trialGridFunctionSpace()
const;
376 const GFSV& testGridFunctionSpace()
const;
377 typename GFSU::Traits::SizeType globalSizeU ()
const;
378 typename GFSV::Traits::SizeType globalSizeV ()
const;
412 template<
typename Gr
idOperatorTuple>
413 static void setupGridOperators(GridOperatorTuple& tuple);
522 #endif // DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH const IG & ig
Definition: constraints.hh:148
void assemble(LocalAssemblerEngine &local_assembler_engine)
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:33
Definition: common/assembler.hh:323
The global assembler which performs the traversing of the integration parts.
Definition: common/assembler.hh:22
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
Definition: common/assembler.hh:325
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, AssemblerInterface, LocalAssemblerInterface > Traits
The traits class.
Definition: common/assembler.hh:351
LocalAssemblerInterface LocalAssembler
The type of the local assembler.
Definition: common/assembler.hh:36
The local assembler which provides the engines that drive the global assembler.
Definition: common/assembler.hh:275
The grid operator represents an operator mapping which corresponds to the (possibly nonlinear) algebr...
Definition: common/assembler.hh:346
Definition: common/assembler.hh:324
The local assembler engine which handles the integration parts as provided by the global assemblers...
Definition: common/assembler.hh:33
Definition: common/assembler.hh:326
const P & p
Definition: constraints.hh:147
Base class for local assembler.
Definition: assemblerutilities.hh:182
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:191