Go to the documentation of this file. 1 #ifndef DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
2 #define DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
4 #include <assemblerutilties.hh>
24 template<
class LocalAssemblerEngine>
89 template<
typename EG,
typename LFSU,
typename LFSV>
95 template<
typename EG,
typename LFSV>
101 template<
typename IG,
typename LFSU_S,
typename LFSV_S,
typename LFSU_N,
typename LFSV_N>
103 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
108 template<
typename IG,
typename LFSV_S,
typename LFSV_N>
114 template<
typename IG,
typename LFSU_S,
typename LFSV_S>
120 template<
typename IG,
typename LFSV_S>
128 template<
typename IG,
typename LFSU_S,
typename LFSV_S>
136 template<
typename IG,
typename 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>
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>
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>
163 template<
typename EG,
typename LFSV>
185 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
187 const LFSV_S & lfsv_s);
189 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
191 const LFSV_S & lfsv_s);
193 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
194 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
196 const LFSV_S & lfsv_s,
197 const LFSV_N & lfsv_n);
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);
203 const LFSV_S & lfsv_s,
204 const LFSV_N & lfsv_n,
205 const LFSV_Coupling & lfsv_coupling);
208 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
210 const LFSV_S & lfsv_s);
212 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s);
214 const LFSV_S & lfsv_s);
216 const LFSU_S & lfsu_s,
const LFSV_S & lfsv_s,
217 const LFSU_N & lfsu_n,
const LFSV_N & lfsv_n);
219 const LFSV_S & lfsv_s,
220 const LFSV_N & lfsv_n);
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);
226 const LFSV_S & lfsv_s,
227 const LFSV_N & lfsv_n,
228 const LFSV_Coupling & lfsv_coupling);
274 template<
typename B,
typename CU,
typename CV>
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>
344 template<
typename GFSU,
typename GFSV,
345 typename MB,
typename DF,
typename RF,
typename JF>
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;
377 typename GFSU::Traits::SizeType
globalSizeU ()
const;
378 typename GFSV::Traits::SizeType
globalSizeV ()
const;
412 template<
typename Gr
idOperatorTuple>
522 #endif // DUNE_PDELAB_GRIDOPERATOR_COMMON_ASSEMBLER_HH
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)
void setJacobian(const J &j)
void preAssembly()
Called directly before assembling.
void onUnbindLFSUVInside(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
bool requireVProcessor() const
Definition: common/assembler.hh:325
Base class for local assembler.
Definition: assemblerutilities.hh:182
void residual(const X &x, R &r) const
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(A &a, const X &x)
LocalAssemblerInterface LocalAssembler
The type of the local assembler.
Definition: common/assembler.hh:36
void assembleUVVolumePostSkeleton(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
void onBindLFSVInside(const IG &ig, const LFSV_S &lfsv_s)
bool requireVSkeleton() const
void onUnbindLFSVInside(const IG &ig, const LFSV_S &lfsv_s)
The global assembler which performs the traversing of the integration parts.
Definition: common/assembler.hh:22
bool requireUVSkeleton() const
void loadCoefficientsLFSUCoupling(const LFSU_Coupling &lfsu_coupling)
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
void loadCoefficientsLFSUOutside(const LFSU_N &lfsu_n)
bool requireSkeleton() const
bool requireUVProcessor() const
LocalResidualAssemblerEngine & localResidualAssemblerEngine(R &r, const X &x)
void setPattern(const P &p)
void postAssembly()
Called last thing after assembling.
void setSolution(const X &x)
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)
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:33
void onUnbindLFSUVCoupling(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n, const LFSU_Coupling &lfsu_coupling, const LFSV_Coupling &lfsv_coupling)
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)
The local assembler which provides the engines that drive the global assembler.
Definition: common/assembler.hh:275
LocalPatternAssemblerEngine & localPatternAssemblerEngine(P &p)
void assembleUVEnrichedCoupling(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n, const LFSU_C &lfsu_c, const LFSV_C &lfsv_c)
void onUnbindLFSVCoupling(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n, const LFSV_Coupling &lfsv_coupling)
void loadCoefficientsLFSUInside(const LFSU_S &lfsu_s)
bool requireUVVolume() const
LocalResidualJacobianAssemblerEngine & localResidualJacobianAssemblerEngine(R &r, A &a, const X &x)
void assembleVEnrichedCoupling(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n, const LFSV_C &lfsv_c)
GFSV::Traits::SizeType globalSizeV() const
bool requireUVVolumePostSkeleton() const
void assembleVVolume(const EG &eg, const LFSV &lfsv)
static void setupGridOperators(GridOperatorTuple &tuple)
const IG & ig
Definition: constraints.hh:148
bool requireVEnrichedCoupling() const
Definition: common/assembler.hh:326
void assembleVVolumePostSkeleton(const EG &eg, const LFSV &lfsv)
GFSU::Traits::SizeType globalSizeU() const
const GFSV & testGridFunctionSpace() const
void onBindLFSUVInside(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
Definition: common/assembler.hh:324
GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, AssemblerInterface, LocalAssemblerInterface > Traits
The traits class.
Definition: common/assembler.hh:351
void interpolate(const typename Traits::Domain &xold, const F &f, typename Traits::Domain &xnew)
Interpolate xnew from f, taking unconstrained values from xold.
The local assembler engine which handles the integration parts as provided by the global assemblers.
Definition: common/assembler.hh:33
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
bool requireUVEnrichedCoupling() const
void fill_pattern(P &globalpattern) const
Determines the sparsity pattern of the jacobian matrix.
void postStep()
Notify local assembler about completion of time step.
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 setTime(TT time)
Set current time of assembling.
void postStage()
Notify local assembler about completion of time step stage.
void assembleUVBoundary(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s)
bool requireVVolume() const
Definition: common/assembler.hh:323
void assembleVProcessor(const IG &ig, const LFSV_S &lfsv_s)
TT suggestTimestep(TT dt) const
Suggest a valid time step size.
void onUnbindLFSV(const EG &eg, const LFSV_S &lfsv_s)
void preStep(TT time, TT dt, std::size_t stages)
Notify local assembler about upcoming time step.
bool requireUVBoundary() const
bool requireVBoundary() const
void setResidual(const R &r)
bool requireVVolumePostSkeleton() const
bool assembleCell(const EG &eg)
void onBindLFSVCoupling(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n, const LFSV_Coupling &lfsv_coupling)
LocalAssemblerInterface & localAssembler()
void assembleUVVolume(const EG &eg, const LFSU &lfsu, const LFSV &lfsv)
const GFSU & trialGridFunctionSpace() const
const LocalAssembler & localAssembler()
Access to the superior local assembler object.
void onBindLFSUVCoupling(const IG &ig, const LFSU_S &lfsu_s, const LFSV_S &lfsv_s, const LFSU_N &lfsu_n, const LFSV_N &lfsv_n const LFSU_Coupling &lfsu_coupling, const LFSV_Coupling &lfsv_coupling)
void assembleVSkeleton(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void preStage(TT time, std::size_t stage)
Notify local assembler about upcoming time step stage.
void assembleVBoundary(const IG &ig, const LFSV_S &lfsv_s)
The grid operator represents an operator mapping which corresponds to the (possibly nonlinear) algebr...
Definition: common/assembler.hh:346
void onBindLFSVOutside(const IG &ig, const LFSV_S &lfsv_s, const LFSV_N &lfsv_n)
void assemble(LocalAssemblerEngine &local_assembler_engine)
void jacobian(const X &x, A &a) const
void setWeight(RF weight)
Set current weight of assembling.
void onBindLFSV(const EG &eg, const LFSV_S &lfsv_s)
const P & p
Definition: constraints.hh:147