DUNE PDELab (git)

assembler.hh
1#ifndef DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
2#define DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
3
4#include <gridoperatorutilities.hh>
5#include <assemblerutilties.hh>
6
7namespace Dune {
8 namespace PDELab {
9
13
24 public:
25 template<class LocalAssemblerEngine>
26 void assemble(LocalAssemblerEngine & local_assembler_engine);
27 };
28
29
35 public:
38
41
49 bool requireSkeleton() const;
50 bool requireSkeletonTwoSided() const;
51 bool requireUVVolume() const;
52 bool requireVVolume() const;
53 bool requireUVSkeleton() const;
54 bool requireVSkeleton() const;
55 bool requireUVBoundary() const;
56 bool requireVBoundary() const;
57 bool requireUVProcessor() const;
58 bool requireVProcessor() const;
59 bool requireUVEnrichedCoupling() const;
60 bool requireVEnrichedCoupling() const;
61 bool requireUVVolumePostSkeleton() const;
62 bool requireVVolumePostSkeleton() const;
64
76 template<typename EG>
77 [[deprecated]] bool assembleCell(const EG & eg)
78 {
79 return false;
80 }
81
91 template<typename EG>
92 bool skipEntity(const EG & eg);
93
103 template<typename IG>
104 bool skipIntersection(const IG & ig);
105
109 template<typename EG, typename LFSU, typename LFSV>
110 void assembleUVVolume(const EG & eg, const LFSU & lfsu, const LFSV & lfsv);
111
115 template<typename EG, typename LFSV>
116 void assembleVVolume(const EG & eg, const LFSV & lfsv);
117
121 template<typename IG, typename LFSU_S, typename LFSV_S, typename LFSU_N, typename LFSV_N>
122 void assembleUVSkeleton(const IG & ig, const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
123 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n);
124
128 template<typename IG, typename LFSV_S, typename LFSV_N>
129 void assembleVSkeleton(const IG & ig, const LFSV_S & lfsv_s, const LFSV_N & lfsv_n);
130
134 template<typename IG, typename LFSU_S, typename LFSV_S>
135 void assembleUVBoundary(const IG & ig, const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
136
140 template<typename IG, typename LFSV_S>
141 void assembleVBoundary(const IG & ig, const LFSV_S & lfsv_s);
142
148 template<typename IG, typename LFSU_S, typename LFSV_S>
149 void assembleUVProcessor(const IG & ig, const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
150
156 template<typename IG, typename LFSV_S>
157 void assembleVProcessor(const IG & ig, const LFSV_S & lfsv_s);
158
159 template<typename IG, typename LFSU_S, typename LFSV_S, typename LFSU_N, typename LFSV_N,
160 typename LFSU_C, typename LFSV_C>
161 void assembleUVEnrichedCoupling(const IG & ig,
162 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
163 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n,
164 const LFSU_C & lfsu_c, const LFSV_C & lfsv_c);
165
166 template<typename IG, typename LFSV_S, typename LFSV_N, typename LFSV_C>
167 void assembleVEnrichedCoupling(const IG & ig,
168 const LFSV_S & lfsv_s,
169 const LFSV_N & lfsv_n,
170 const LFSV_C & lfsv_c);
171
176 template<typename EG, typename LFSU, typename LFSV>
177 void assembleUVVolumePostSkeleton(const EG & eg, const LFSU & lfsu, const LFSV & lfsv);
178
183 template<typename EG, typename LFSV>
184 void assembleVVolumePostSkeleton(const EG & eg, const LFSV & lfsv);
185
188
191
204 void onBindLFSUV(const EG & eg,
205 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
206 void onBindLFSV(const EG & eg,
207 const LFSV_S & lfsv_s);
208 void onBindLFSUVInside(const IG & ig,
209 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
210 void onBindLFSVInside(const IG & ig,
211 const LFSV_S & lfsv_s);
212 void onBindLFSUVOutside(const IG & ig,
213 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
214 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n);
215 void onBindLFSVOutside(const IG & ig,
216 const LFSV_S & lfsv_s,
217 const LFSV_N & lfsv_n);
218 void onBindLFSUVCoupling(const IG & ig,
219 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
220 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n
221 const LFSU_Coupling & lfsu_coupling, const LFSV_Coupling & lfsv_coupling);
222 void onBindLFSVCoupling(const IG & ig,
223 const LFSV_S & lfsv_s,
224 const LFSV_N & lfsv_n,
225 const LFSV_Coupling & lfsv_coupling);
226
227 void onUnbindLFSUV(const EG & eg,
228 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
229 void onUnbindLFSV(const EG & eg,
230 const LFSV_S & lfsv_s);
231 void onUnbindLFSUVInside(const IG & ig,
232 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s);
233 void onUnbindLFSVInside(const IG & ig,
234 const LFSV_S & lfsv_s);
235 void onUnbindLFSUVOutside(const IG & ig,
236 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
237 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n);
238 void onUnbindLFSVOutside(const IG & ig,
239 const LFSV_S & lfsv_s,
240 const LFSV_N & lfsv_n);
241 void onUnbindLFSUVCoupling(const IG & ig,
242 const LFSU_S & lfsu_s, const LFSV_S & lfsv_s,
243 const LFSU_N & lfsu_n, const LFSV_N & lfsv_n,
244 const LFSU_Coupling & lfsu_coupling, const LFSV_Coupling & lfsv_coupling);
245 void onUnbindLFSVCoupling(const IG & ig,
246 const LFSV_S & lfsv_s,
247 const LFSV_N & lfsv_n,
248 const LFSV_Coupling & lfsv_coupling);
249
260 void loadCoefficientsLFSUInside(const LFSU_S & lfsu_s);
261 void loadCoefficientsLFSUOutside(const LFSU_N & lfsu_n);
262 void loadCoefficientsLFSUCoupling(const LFSU_Coupling & lfsu_coupling);
275 void setSolution(const X& x);
276 void setPattern(const P& p);
277 void setJacobian(const J & j);
278 void setResidual(const R& r);
281 };
282
294 template<typename B, typename CU, typename CV>
296 {
297 public:
298
304 template<class TT>
305 void setTime(TT time);
306
308 template<typename TT>
309 void preStep (TT time, TT dt, std::size_t stages);
310
312 void postStep ();
313
315 template<typename TT>
316 void preStage (TT time, std::size_t stage);
317
319 void postStage ();
320
322 template<typename TT>
323 TT suggestTimestep (TT dt) const;
324
328 template<class RF>
329 void setWeight(RF weight);
330
334 LocalPatternAssemblerEngine & localPatternAssemblerEngine(P & p);
335 LocalResidualAssemblerEngine & localResidualAssemblerEngine(R & r, const X & x);
336 LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(A & a, const X & x);
337 LocalResidualJacobianAssemblerEngine & localResidualJacobianAssemblerEngine(R & r, A & a, const X & x);
343 class LocalPatternAssemblerEngine : public LocalAssemblerEngine {};
344 class LocalResidualAssemblerEngine : public LocalAssemblerEngine {};
345 class LocalJacobianAssemblerEngine : public LocalAssemblerEngine {};
346 class LocalResidualJacobianAssemblerEngine : public LocalAssemblerEngine {};
349 };
350
364 template<typename GFSU, typename GFSV,
365 typename MB, typename DF, typename RF, typename JF>
367 public:
368
370 typedef GridOperatorTraits
371 <GFSU,GFSV,MB,DF,RF,JF,CU,CV,AssemblerInterface,LocalAssemblerInterface> Traits;
372
374 template<typename P>
375 void fill_pattern (P& globalpattern) const;
376
379 template<typename X, typename R>
380 void residual (const X& x, R& r) const;
381
384 template<typename X, typename A>
385 void jacobian (const X& x, A& a) const;
386
389 Assembler & assembler();
390 LocalAssemblerInterface & localAssembler();
392
395 const GFSU& trialGridFunctionSpace() const;
396 const GFSV& testGridFunctionSpace() const;
397 typename GFSU::Traits::SizeType globalSizeU () const;
398 typename GFSV::Traits::SizeType globalSizeV () const;
400
402
408 template<typename F>
409 void interpolate(const typename Traits::Domain& xold,
410 const F& f,
411 typename Traits::Domain& xnew);
412
417
432 template<typename GridOperatorTuple>
433 static void setupGridOperators(GridOperatorTuple& tuple);
434
435 };
437 };
438};
439
542#endif // DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
The global assembler which performs the traversing of the integration parts.
Definition: assembler.hh:23
The grid operator represents an operator mapping which corresponds to the (possibly nonlinear) algebr...
Definition: assembler.hh:366
GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, AssemblerInterface, LocalAssemblerInterface > Traits
The traits class.
Definition: assembler.hh:371
void interpolate(const typename Traits::Domain &xold, const F &f, typename Traits::Domain &xnew)
Interpolate xnew from f, taking unconstrained values from xold.
void residual(const X &x, R &r) const
static void setupGridOperators(GridOperatorTuple &tuple)
void jacobian(const X &x, A &a) const
void fill_pattern(P &globalpattern) const
Determines the sparsity pattern of the jacobian matrix.
Base class for local assembler.
Definition: assemblerutilities.hh:217
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)
const LocalAssembler & localAssembler()
Access to the superior local assembler object.
void preAssembly()
Called directly before assembling.
void assembleVVolumePostSkeleton(const EG &eg, const LFSV &lfsv)
bool assembleCell(const EG &eg)
Deprecated. Use skipEntity instead.
Definition: assembler.hh:77
void assembleUVVolumePostSkeleton(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
LocalAssemblerInterface LocalAssembler
The type of the local assembler.
Definition: assembler.hh:37
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)
The local assembler which provides the engines that drive the global assembler.
Definition: assembler.hh:296
void preStep(TT time, TT dt, std::size_t stages)
Notify local assembler about upcoming time step.
void setWeight(RF weight)
Set current weight of assembling.
void postStage()
Notify local assembler about completion of time step stage.
void postStep()
Notify local assembler about completion of time step.
void preStage(TT time, std::size_t stage)
Notify local assembler about upcoming time step stage.
void setTime(TT time)
Set current time of assembling.
TT suggestTimestep(TT dt) const
Suggest a valid time step size.
Dune namespace.
Definition: alignedallocator.hh:13
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:34
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)