3 #ifndef DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH 4 #define DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH 16 template<
typename GO0,
typename GO1,
bool implicit = true>
41 <
typename GO0::Traits::TrialGridFunctionSpace,
42 typename GO0::Traits::TestGridFunctionSpace,
43 typename GO0::Traits::MatrixBackend,
44 typename GO0::Traits::DomainField,
45 typename GO0::Traits::RangeField,
46 typename GO0::Traits::JacobianField,
47 typename GO0::Traits::TrialGridFunctionSpaceConstraints,
48 typename GO0::Traits::TestGridFunctionSpaceConstraints,
59 template <
typename MFT>
77 local_assembler(la0,la1, const_residual)
79 GO0::setupGridOperators(std::tie(go0_,go1_));
90 DUNE_THROW(Dune::Exception,
"This function should not be called in explicit mode");
96 DUNE_THROW(Dune::Exception,
"This function should not be called in explicit mode");
103 return global_assembler.trialGridFunctionSpace();
109 return global_assembler.testGridFunctionSpace();
113 typename Traits::TrialGridFunctionSpace::Traits::SizeType
globalSizeU ()
const 119 typename Traits::TestGridFunctionSpace::Traits::SizeType
globalSizeV ()
const 124 Assembler &
assembler()
const {
return global_assembler; }
134 PatternEngine & pattern_engine = local_assembler.localPatternAssemblerEngine(p);
135 global_assembler.assemble(pattern_engine);
140 PatternEngine & pattern_engine = local_assembler.localExplicitPatternAssemblerEngine(p);
141 global_assembler.assemble(pattern_engine);
146 void preStage(
unsigned int stage,
const std::vector<Domain*> & x)
149 DUNE_THROW(Dune::Exception,
"This function should not be called in explicit mode");
152 local_assembler.setStage(stage);
153 PreStageEngine & prestage_engine = local_assembler.localPreStageAssemblerEngine(x);
154 global_assembler.assemble(prestage_engine);
163 DUNE_THROW(Dune::Exception,
"This function should not be called in explicit mode");
166 ResidualEngine & residual_engine = local_assembler.localResidualAssemblerEngine(r,x);
167 global_assembler.assemble(residual_engine);
173 void jacobian(
const Domain & x, Jacobian & a)
const 176 DUNE_THROW(Dune::Exception,
"This function should not be called in explicit mode");
179 JacobianEngine & jacobian_engine = local_assembler.localJacobianAssemblerEngine(a,x);
180 global_assembler.assemble(jacobian_engine);
187 Jacobian & a, Range & r1, Range & r0)
190 DUNE_THROW(Dune::Exception,
"This function should not be called in implicit mode");
192 local_assembler.setStage(stage);
195 ExplicitJacobianResidualEngine;
197 ExplicitJacobianResidualEngine & jacobian_residual_engine
198 = local_assembler.localExplicitJacobianResidualAssemblerEngine(a,r0,r1,x);
200 global_assembler.assemble(jacobian_residual_engine);
204 template<
typename F,
typename X>
205 void interpolate (
unsigned stage,
const X& xold, F& f, X& x)
const 208 f.setTime(local_assembler.timeAtStage(stage));
210 go0.localAssembler().setTime(local_assembler.timeAtStage(stage));
213 go0.interpolate(xold,f,x);
222 local_assembler.setMethod(method_);
228 local_assembler.setMethod(method_);
229 local_assembler.preStep(time_,dt_,method_.
s());
249 Real suggested_dt = std::min(la0.suggestTimestep(dt),la1.suggestTimestep(dt));
259 const_residual =
Range(go0.testGridFunctionSpace());
264 go0.make_consistent(a);
269 return go0.matrixBackend();
274 return local_assembler.trialConstraints();
278 Assembler & global_assembler;
281 LocalAssemblerDT0 & la0;
282 LocalAssemblerDT1 & la1;
283 Range const_residual;
284 mutable LocalAssembler local_assembler;
289 #endif // DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH
LocalAssembler::Real Real
The type for real number e.g. time.
Definition: gridoperator/onestep.hh:66
void explicit_jacobian_residual(unsigned int stage, const std::vector< Domain *> &x, Jacobian &a, Range &r1, Range &r0)
Assemble jacobian and residual simultaneously for explicit treatment.
Definition: gridoperator/onestep.hh:186
GFSU TrialGridFunctionSpace
The trial grid function space.
Definition: gridoperatorutilities.hh:37
Assembler & assembler() const
Definition: gridoperator/onestep.hh:124
Traits::TrialGridFunctionSpace::Traits::SizeType globalSizeU() const
Get dimension of space u.
Definition: gridoperator/onestep.hh:113
void preStage(unsigned int stage, const std::vector< Domain *> &x)
Assemble constant part of residual.
Definition: gridoperator/onestep.hh:146
Definition: gridoperator/onestep.hh:17
GO::Traits::TrialGridFunctionSpaceConstraints TrialGridFunctionSpaceConstraints
The type of the trial grid function space constraints.
Definition: assemblerutilities.hh:33
void postStep()
to be called after step is completed
Definition: gridoperator/onestep.hh:233
Dune::PDELab::GridOperatorTraits< typename GO0::Traits::TrialGridFunctionSpace, typename GO0::Traits::TestGridFunctionSpace, typename GO0::Traits::MatrixBackend, typename GO0::Traits::DomainField, typename GO0::Traits::RangeField, typename GO0::Traits::JacobianField, typename GO0::Traits::TrialGridFunctionSpaceConstraints, typename GO0::Traits::TestGridFunctionSpaceConstraints, Assembler, LocalAssembler > Traits
The grid operator traits.
Definition: gridoperator/onestep.hh:50
GO0::Traits::Assembler Assembler
The global UDG assembler type.
Definition: gridoperator/onestep.hh:25
void fill_pattern(Pattern &p) const
Fill pattern of jacobian matrix.
Definition: gridoperator/onestep.hh:129
Jacobian Type
Definition: gridoperator/onestep.hh:62
OneStepLocalAssembler< OneStepGridOperator, LocalAssemblerDT0, LocalAssemblerDT1 > LocalAssembler
The local assembler type.
Definition: gridoperator/onestep.hh:34
void interpolate(unsigned stage, const X &xold, F &f, X &x) const
Interpolate constrained values from given function f.
Definition: gridoperator/onestep.hh:205
void preStep(const TimeSteppingParameterInterface< Real > &method_, Real time_, Real dt_)
parametrize assembler with a time-stepping method
Definition: gridoperator/onestep.hh:226
Definition: onestep/localassembler.hh:147
OneStepGridOperator(GO0 &go0_, GO1 &go1_)
Constructor for non trivial constraints.
Definition: gridoperator/onestep.hh:72
Base parameter class for time stepping scheme parameters.
Definition: onestepparameter.hh:23
void update()
Definition: gridoperator/onestep.hh:255
void setMethod(const TimeSteppingParameterInterface< Real > &method_)
set time stepping method
Definition: gridoperator/onestep.hh:220
Dune::PDELab::Backend::Vector< GFSV, RF > Range
The type of the range (residual).
Definition: gridoperatorutilities.hh:65
GO0::Traits::LocalAssembler LocalAssemblerDT0
Definition: gridoperator/onestep.hh:29
GO0::Pattern Pattern
The sparsity pattern container for the jacobian matrix.
Definition: gridoperator/onestep.hh:22
void jacobian(const Domain &x, Jacobian &a) const
Assemble jacobian.
Definition: gridoperator/onestep.hh:173
MB MatrixBackend
The matrix backend of the grid operator.
Definition: gridoperatorutilities.hh:51
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints trialConstraints() const
Definition: gridoperator/onestep.hh:272
void residual(const Domain &x, Range &r) const
Assemble residual.
Definition: gridoperator/onestep.hh:160
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Definition: onestep/localassembler.hh:147
LocalAssemblerDT1 ::LocalPatternAssemblerEngine LocalExplicitPatternAssemblerEngine
Definition: onestep/localassembler.hh:55
Traits::TestGridFunctionSpace::Traits::SizeType globalSizeV() const
Get dimension of space v.
Definition: gridoperator/onestep.hh:119
GO0::BorderDOFExchanger BorderDOFExchanger
The BorderDOFExchanger.
Definition: gridoperator/onestep.hh:37
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
LocalAssembler & localAssembler() const
Definition: gridoperator/onestep.hh:126
void divideMassTermByDeltaT()
Definition: gridoperator/onestep.hh:87
GFSV TestGridFunctionSpace
The test grid function space.
Definition: gridoperatorutilities.hh:40
Traits::Range Range
Definition: gridoperator/onestep.hh:55
void multiplySpatialTermByDeltaT()
Definition: gridoperator/onestep.hh:93
Traits::Jacobian Jacobian
Definition: gridoperator/onestep.hh:56
GO1::Traits::LocalAssembler LocalAssemblerDT1
Definition: gridoperator/onestep.hh:30
Definition: onestep/localassembler.hh:147
LocalAssembler::OneStepParameters OneStepParameters
The type of the one step method parameters.
Definition: gridoperator/onestep.hh:69
const Traits::TrialGridFunctionSpace & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: gridoperator/onestep.hh:101
virtual unsigned s() const =0
Return number of stages of the method.
const Traits::TestGridFunctionSpace & testGridFunctionSpace() const
Get the test grid function space.
Definition: gridoperator/onestep.hh:107
const Traits::MatrixBackend & matrixBackend() const
Definition: gridoperator/onestep.hh:267
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:989
Real suggestTimestep(Real dt) const
to be called once before each stage
Definition: gridoperator/onestep.hh:247
const P & p
Definition: constraints.hh:147
void postStage()
to be called after stage is completed
Definition: gridoperator/onestep.hh:240
void make_consistent(Jacobian &a) const
Definition: gridoperator/onestep.hh:262
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:33
Traits::Domain Domain
Definition: gridoperator/onestep.hh:54
Definition: gridoperator/onestep.hh:60
Traits::RangeField Real
The local operators type for real numbers e.g. time.
Definition: onestep/localassembler.hh:87