dune-pdelab 2.7-git
Loading...
Searching...
No Matches
gridoperator/onestep.hh
Go to the documentation of this file.
1// -*- tab-width: 2; indent-tabs-mode: nil -*-
2// vi: set et ts=2 sw=2 sts=2:
3#ifndef DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH
4#define DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH
5
6#include <tuple>
7
13
14namespace Dune{
15 namespace PDELab{
16
17 template<typename GO0, typename GO1, bool implicit = true>
19 {
20 public:
21
23 typedef typename GO0::Pattern Pattern;
24
26 typedef typename GO0::Traits::Assembler Assembler;
27
30 typedef typename GO0::Traits::LocalAssembler LocalAssemblerDT0;
31 typedef typename GO1::Traits::LocalAssembler LocalAssemblerDT1;
33
36
38 typedef typename GO0::BorderDOFExchanger BorderDOFExchanger;
39
42 <typename GO0::Traits::TrialGridFunctionSpace,
43 typename GO0::Traits::TestGridFunctionSpace,
44 typename GO0::Traits::MatrixBackend,
45 typename GO0::Traits::DomainField,
46 typename GO0::Traits::RangeField,
47 typename GO0::Traits::JacobianField,
48 typename GO0::Traits::TrialGridFunctionSpaceConstraints,
49 typename GO0::Traits::TestGridFunctionSpaceConstraints,
52
55 typedef typename Traits::Domain Domain;
56 typedef typename Traits::Range Range;
57 typedef typename Traits::Jacobian Jacobian;
59
60 template <typename MFT>
62 {
63 typedef Jacobian Type;
64 };
65
67 typedef typename LocalAssembler::Real Real;
68
71
73 OneStepGridOperator(GO0 & go0_, GO1 & go1_)
74 : global_assembler(go0_.assembler()),
75 go0(go0_), go1(go1_),
76 la0(go0_.localAssembler()), la1(go1_.localAssembler()),
77 const_residual( go0_.testGridFunctionSpace() ),
78 local_assembler(la0,la1, const_residual)
79 {
80 GO0::setupGridOperators(std::tie(go0_,go1_));
81 if(not implicit)
83 }
84
89 {
90 if(not implicit)
91 DUNE_THROW(Dune::Exception,"This function should not be called in explicit mode");
93 }
95 {
96 if(not implicit)
97 DUNE_THROW(Dune::Exception,"This function should not be called in explicit mode");
99 }
100
103 {
104 return global_assembler.trialGridFunctionSpace();
105 }
106
109 {
110 return global_assembler.testGridFunctionSpace();
111 }
112
114 typename Traits::TrialGridFunctionSpace::Traits::SizeType globalSizeU () const
115 {
116 return trialGridFunctionSpace().globalSize();
117 }
118
120 typename Traits::TestGridFunctionSpace::Traits::SizeType globalSizeV () const
121 {
122 return testGridFunctionSpace().globalSize();
123 }
124
125 Assembler & assembler() const { return global_assembler; }
126
127 LocalAssembler & localAssembler() const { return local_assembler; }
128
130 void fill_pattern(Pattern & p) const
131 {
132 if(implicit)
133 {
134 typedef typename LocalAssembler::LocalPatternAssemblerEngine PatternEngine;
135 PatternEngine & pattern_engine = local_assembler.localPatternAssemblerEngine(p);
136 global_assembler.assemble(pattern_engine);
137 }
138 else
139 {
140 typedef typename LocalAssembler::LocalExplicitPatternAssemblerEngine PatternEngine;
141 PatternEngine & pattern_engine = local_assembler.localExplicitPatternAssemblerEngine(p);
142 global_assembler.assemble(pattern_engine);
143 }
144 }
145
147 void preStage(unsigned int stage, const std::vector<Domain*> & x)
148 {
149 if(not implicit)
150 DUNE_THROW(Dune::Exception,"This function should not be called in explicit mode");
151
152 typedef typename LocalAssembler::LocalPreStageAssemblerEngine PreStageEngine;
153 local_assembler.setStage(stage);
154 PreStageEngine & prestage_engine = local_assembler.localPreStageAssemblerEngine(x);
155 global_assembler.assemble(prestage_engine);
156 }
157
159 void residual(const Domain & x, Range & r) const
160 {
161 if(not implicit)
162 DUNE_THROW(Dune::Exception,"This function should not be called in explicit mode");
163
164 typedef typename LocalAssembler::LocalResidualAssemblerEngine ResidualEngine;
165 ResidualEngine & residual_engine = local_assembler.localResidualAssemblerEngine(r,x);
166 global_assembler.assemble(residual_engine);
167 }
168
170 void jacobian(const Domain & x, Jacobian & a) const
171 {
172 if(not implicit)
173 DUNE_THROW(Dune::Exception,"This function should not be called in explicit mode");
174
175 typedef typename LocalAssembler::LocalJacobianAssemblerEngine JacobianEngine;
176 JacobianEngine & jacobian_engine = local_assembler.localJacobianAssemblerEngine(a,x);
177 global_assembler.assemble(jacobian_engine);
178 }
179
181 void explicit_jacobian_residual(unsigned int stage, const std::vector<Domain*> & x,
182 Jacobian & a, Range & r1, Range & r0)
183 {
184 if(implicit)
185 DUNE_THROW(Dune::Exception,"This function should not be called in implicit mode");
186
187 local_assembler.setStage(stage);
188
190 ExplicitJacobianResidualEngine;
191
192 ExplicitJacobianResidualEngine & jacobian_residual_engine
193 = local_assembler.localExplicitJacobianResidualAssemblerEngine(a,r0,r1,x);
194 global_assembler.assemble(jacobian_residual_engine);
195 }
196
198 void jacobian_apply(const Domain & update, Range & result) const
199 {
200 if ((not la0.localOperator().isLinear) or (not la1.localOperator().isLinear))
201 DUNE_THROW(Dune::Exception, "Your trying to use a linear jacobian apply for a non linear problem.");
202 global_assembler.assemble(local_assembler.localJacobianApplyAssemblerEngine(update, result));
203 }
204
206 void jacobian_apply(const Domain & solution, const Domain & update, Range & result) const
207 {
208 if (la0.localOperator().isLinear and la1.localOperator().isLinear)
209 DUNE_THROW(Dune::Exception, "Your trying to use a non linear jacobian apply for a linear problem.");
210 global_assembler.assemble(local_assembler.localJacobianApplyAssemblerEngine(solution, update, result));
211 }
212
214 template<typename F, typename X>
215 void interpolate (unsigned stage, const X& xold, F& f, X& x) const
216 {
217 // Set time in boundary value function
218 f.setTime(local_assembler.timeAtStage(stage));
219
220 go0.localAssembler().setTime(local_assembler.timeAtStage(stage));
221
222 // Interpolate
223 go0.interpolate(xold,f,x);
224
225 // Copy non-constrained dofs from old time step
227
228 // If there are Dirichlet-constrained DOFs on processor boundaries, we need to make
229 // those DOFs consistent again
230 auto es = trialGridFunctionSpace().entitySet();
231 if (es.comm().size() > 1)
232 {
234 es.communicate(data_handle,InteriorBorder_All_Interface,ForwardCommunication);
235 }
236 }
237
240 {
241 local_assembler.setMethod(method_);
242 }
243
245 void preStep (const TimeSteppingParameterInterface<Real>& method_, Real time_, Real dt_)
246 {
247 local_assembler.setMethod(method_);
248 local_assembler.preStep(time_,dt_,method_.s());
249 }
250
252 void postStep ()
253 {
254 la0.postStep();
255 la1.postStep();
256 }
257
259 void postStage ()
260 {
261 la0.postStage();
262 la1.postStage();
263 }
264
267 {
268 Real suggested_dt = std::min(la0.suggestTimestep(dt),la1.suggestTimestep(dt));
269 if (trialGridFunctionSpace().gridView().comm().size()>1)
270 suggested_dt = trialGridFunctionSpace().gridView().comm().min(suggested_dt);
271 return suggested_dt;
272 }
273
274 void update()
275 {
276 go0.update();
277 go1.update();
278 const_residual = Range(go0.testGridFunctionSpace());
279 }
280
282 {
283 go0.make_consistent(a);
284 }
285
286 const typename Traits::MatrixBackend& matrixBackend() const
287 {
288 return go0.matrixBackend();
289 }
290
291 const typename LocalAssembler::Traits::TrialGridFunctionSpaceConstraints trialConstraints() const
292 {
293 return local_assembler.trialConstraints();
294 }
295
296 private:
297 Assembler & global_assembler;
298 GO0 & go0;
299 GO1 & go1;
300 LocalAssemblerDT0 & la0;
301 LocalAssemblerDT1 & la1;
302 Range const_residual;
303 mutable LocalAssembler local_assembler;
304 };
305
306 }
307}
308#endif // DUNE_PDELAB_GRIDOPERATOR_ONESTEP_HH
const P & p
Definition: constraints.hh:148
void copy_nonconstrained_dofs(const CG &cg, const XG &xgin, XG &xgout)
Definition: constraints.hh:987
virtual unsigned s() const =0
Return number of stages of the method.
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Definition: genericdatahandle.hh:730
const CU & trialConstraints() const
get the constraints on the trial grid function space
Definition: assemblerutilities.hh:233
Traits class for the grid operator.
Definition: gridoperatorutilities.hh:34
GFSU TrialGridFunctionSpace
The trial grid function space.
Definition: gridoperatorutilities.hh:37
GFSV TestGridFunctionSpace
The test grid function space.
Definition: gridoperatorutilities.hh:40
MB MatrixBackend
The matrix backend of the grid operator.
Definition: gridoperatorutilities.hh:51
Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF > Jacobian
The type of the jacobian.
Definition: gridoperatorutilities.hh:72
Dune::PDELab::Backend::Vector< GFSV, RF > Range
The type of the range (residual).
Definition: gridoperatorutilities.hh:65
Dune::PDELab::Backend::Vector< GFSU, DF > Domain
The type of the domain (solution).
Definition: gridoperatorutilities.hh:58
Definition: gridoperator/onestep.hh:19
Traits::Domain Domain
Definition: gridoperator/onestep.hh:55
void postStep()
to be called after step is completed
Definition: gridoperator/onestep.hh:252
const LocalAssembler::Traits::TrialGridFunctionSpaceConstraints trialConstraints() const
Definition: gridoperator/onestep.hh:291
void interpolate(unsigned stage, const X &xold, F &f, X &x) const
Interpolate constrained values from given function f.
Definition: gridoperator/onestep.hh:215
const Traits::TestGridFunctionSpace & testGridFunctionSpace() const
Get the test grid function space.
Definition: gridoperator/onestep.hh:108
void divideMassTermByDeltaT()
Definition: gridoperator/onestep.hh:88
LocalAssembler & localAssembler() const
Definition: gridoperator/onestep.hh:127
void preStage(unsigned int stage, const std::vector< Domain * > &x)
Assemble constant part of residual.
Definition: gridoperator/onestep.hh:147
GO1::Traits::LocalAssembler LocalAssemblerDT1
Definition: gridoperator/onestep.hh:31
LocalAssembler::Real Real
The type for real number e.g. time.
Definition: gridoperator/onestep.hh:67
void setMethod(const TimeSteppingParameterInterface< Real > &method_)
set time stepping method
Definition: gridoperator/onestep.hh:239
void residual(const Domain &x, Range &r) const
Assemble residual.
Definition: gridoperator/onestep.hh:159
void multiplySpatialTermByDeltaT()
Definition: gridoperator/onestep.hh:94
void fill_pattern(Pattern &p) const
Fill pattern of jacobian matrix.
Definition: gridoperator/onestep.hh:130
Real suggestTimestep(Real dt) const
to be called once before each stage
Definition: gridoperator/onestep.hh:266
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:181
GO0::Traits::LocalAssembler LocalAssemblerDT0
Definition: gridoperator/onestep.hh:30
void jacobian(const Domain &x, Jacobian &a) const
Assemble jacobian.
Definition: gridoperator/onestep.hh:170
void make_consistent(Jacobian &a) const
Definition: gridoperator/onestep.hh:281
OneStepLocalAssembler< OneStepGridOperator, LocalAssemblerDT0, LocalAssemblerDT1 > LocalAssembler
The local assembler type.
Definition: gridoperator/onestep.hh:35
GO0::BorderDOFExchanger BorderDOFExchanger
The BorderDOFExchanger.
Definition: gridoperator/onestep.hh:38
void update()
Definition: gridoperator/onestep.hh:274
void preStep(const TimeSteppingParameterInterface< Real > &method_, Real time_, Real dt_)
parametrize assembler with a time-stepping method
Definition: gridoperator/onestep.hh:245
const Traits::TrialGridFunctionSpace & trialGridFunctionSpace() const
Get the trial grid function space.
Definition: gridoperator/onestep.hh:102
void postStage()
to be called after stage is completed
Definition: gridoperator/onestep.hh:259
Traits::TestGridFunctionSpace::Traits::SizeType globalSizeV() const
Get dimension of space v.
Definition: gridoperator/onestep.hh:120
OneStepGridOperator(GO0 &go0_, GO1 &go1_)
Constructor for non trivial constraints.
Definition: gridoperator/onestep.hh:73
Traits::TrialGridFunctionSpace::Traits::SizeType globalSizeU() const
Get dimension of space u.
Definition: gridoperator/onestep.hh:114
Traits::Range Range
Definition: gridoperator/onestep.hh:56
GO0::Traits::Assembler Assembler
The global UDG assembler type.
Definition: gridoperator/onestep.hh:26
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:51
Assembler & assembler() const
Definition: gridoperator/onestep.hh:125
LocalAssembler::OneStepParameters OneStepParameters
The type of the one step method parameters.
Definition: gridoperator/onestep.hh:70
const Traits::MatrixBackend & matrixBackend() const
Definition: gridoperator/onestep.hh:286
Traits::Jacobian Jacobian
Definition: gridoperator/onestep.hh:57
void jacobian_apply(const Domain &update, Range &result) const
Apply jacobian matrix to the vector update without explicitly assembling it.
Definition: gridoperator/onestep.hh:198
void jacobian_apply(const Domain &solution, const Domain &update, Range &result) const
Apply jacobian matrix to the vector update without explicitly assembling it.
Definition: gridoperator/onestep.hh:206
GO0::Pattern Pattern
The sparsity pattern container for the jacobian matrix.
Definition: gridoperator/onestep.hh:23
Definition: gridoperator/onestep.hh:62
Jacobian Type
Definition: gridoperator/onestep.hh:63
The local assembler for one step methods.
Definition: onestep/localassembler.hh:34
void setStage(int stage_)
Set the current stage of the one step scheme.
Definition: onestep/localassembler.hh:146
LocalPatternAssemblerEngine & localPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: onestep/localassembler.hh:185
LocalPreStageAssemblerEngine & localPreStageAssemblerEngine(const std::vector< typename Traits::Solution * > &x)
Definition: onestep/localassembler.hh:194
Real timeAtStage(int stage_) const
Access time at given stage.
Definition: onestep/localassembler.hh:162
void setMethod(const OneStepParameters &method_)
Set the one step method parameters.
Definition: onestep/localassembler.hh:140
void setDTAssemblingMode(DTAssemblingMode dt_mode_)
Definition: onestep/localassembler.hh:156
LocalResidualAssemblerEngine & localResidualAssemblerEngine(typename Traits::Residual &r, const typename Traits::Solution &x)
Definition: onestep/localassembler.hh:204
LocalJacobianApplyAssemblerEngine & localJacobianApplyAssemblerEngine(const typename Traits::Domain &update, typename Traits::Range &result)
Definition: onestep/localassembler.hh:252
void preStep(Real time_, Real dt_, int stages_)
Definition: onestep/localassembler.hh:109
LA1::LocalPatternAssemblerEngine LocalExplicitPatternAssemblerEngine
Definition: onestep/localassembler.hh:57
LocalExplicitPatternAssemblerEngine & localExplicitPatternAssemblerEngine(typename Traits::MatrixPattern &p)
Definition: onestep/localassembler.hh:225
LocalExplicitJacobianResidualAssemblerEngine & localExplicitJacobianResidualAssemblerEngine(typename Traits::Jacobian &a, typename Traits::Residual &r0, typename Traits::Residual &r1, const std::vector< typename Traits::Solution * > &x)
Definition: onestep/localassembler.hh:233
LocalJacobianAssemblerEngine & localJacobianAssemblerEngine(typename Traits::Jacobian &a, const typename Traits::Solution &x)
Definition: onestep/localassembler.hh:215
Traits::RangeField Real
The local operators type for real numbers e.g. time.
Definition: onestep/localassembler.hh:90
Base parameter class for time stepping scheme parameters.
Definition: onestepparameter.hh:44