version 3.9
multistagefvassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_EXPERIMENTAL_MULTISTAGE_FV_ASSEMBLER_HH
14#define DUMUX_EXPERIMENTAL_MULTISTAGE_FV_ASSEMBLER_HH
15
16#include <vector>
17#include <deque>
18#include <type_traits>
19#include <memory>
20
21#include <dune/istl/matrixindexset.hh>
22
26
30
34
37
41
44
46
47template<class DiscretizationMethod>
49
50template<class DM>
51struct LocalAssemblerChooser<DiscretizationMethods::CVFE<DM>>
52{
53 template<class TypeTag, class Impl, DiffMethod diffMethod>
55};
56
57template<>
58struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa>
59{
60 template<class TypeTag, class Impl, DiffMethod diffMethod>
62};
63
64template<>
65struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa>
66{
67 template<class TypeTag, class Impl, DiffMethod diffMethod>
69};
70
71template<class TypeTag, class Impl, DiffMethod diffMethod>
74>::template type<TypeTag, Impl, diffMethod>;
75
76} // end namespace Dumux::Detail
77
78namespace Dumux::Experimental {
79
87template<class TypeTag, DiffMethod diffMethod>
89{
91 using GridView = typename GridGeo::GridView;
93 using Element = typename GridView::template Codim<0>::Entity;
94 using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed;
95
96 static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box;
97
100
101public:
107
109
110 using GridGeometry = GridGeo;
112
118 MultiStageFVAssembler(std::shared_ptr<const Problem> problem,
119 std::shared_ptr<const GridGeometry> gridGeometry,
120 std::shared_ptr<GridVariables> gridVariables,
121 std::shared_ptr<const Experimental::MultiStageMethod<Scalar>> msMethod,
122 const SolutionVector& prevSol)
123 : timeSteppingMethod_(msMethod)
124 , problem_(problem)
125 , gridGeometry_(gridGeometry)
126 , gridVariables_(gridVariables)
127 , prevSol_(&prevSol)
128 {
130 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
132 && getParam<bool>("Assembly.Multithreading", true);
133
134 maybeComputeColors_();
135 }
136
141 template<class PartialReassembler = DefaultPartialReassembler>
142 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
143 {
144 resetJacobian_(partialReassembler);
145
146 resetResidual_();
147 spatialOperatorEvaluations_.back() = 0.0;
148 temporalOperatorEvaluations_.back() = 0.0;
149
150 if (stageParams_->size() != spatialOperatorEvaluations_.size())
151 DUNE_THROW(Dune::InvalidStateException, "Wrong number of residuals");
152
153 assemble_([&](const Element& element)
154 {
155 LocalAssembler localAssembler(*this, element, curSol);
156 localAssembler.assembleJacobianAndResidual(
157 *jacobian_, *residual_, *gridVariables_,
158 *stageParams_,
159 temporalOperatorEvaluations_.back(),
160 spatialOperatorEvaluations_.back(),
161 constrainedDofs_,
162 partialReassembler
163 );
164 });
165
166 // assemble the full residual for the time integration stage
167 auto constantResidualComponent = (*residual_);
168 constantResidualComponent = 0.0;
169 for (std::size_t k = 0; k < stageParams_->size()-1; ++k)
170 {
171 if (!stageParams_->skipTemporal(k))
172 constantResidualComponent.axpy(stageParams_->temporalWeight(k), temporalOperatorEvaluations_[k]);
173 if (!stageParams_->skipSpatial(k))
174 constantResidualComponent.axpy(stageParams_->spatialWeight(k), spatialOperatorEvaluations_[k]);
175 }
176
177 // masked summation of constant residual component onto this stage's resiudal component
178 for (std::size_t i = 0; i < constantResidualComponent.size(); ++i)
179 for (std::size_t ii = 0; ii < constantResidualComponent[i].size(); ++ii)
180 (*residual_)[i][ii] += constrainedDofs_[i][ii] > 0.5 ? 0.0 : constantResidualComponent[i][ii];
181 }
182
185 { DUNE_THROW(Dune::NotImplemented, "residual"); }
186
192 {
193 jacobian_ = std::make_shared<JacobianMatrix>();
194 jacobian_->setBuildMode(JacobianMatrix::random);
195 residual_ = std::make_shared<ResidualType>();
196
197 setResidualSize_(*residual_);
198 setJacobianPattern_();
199 }
200
205 {
206 setResidualSize_(*residual_);
207 setJacobianPattern_();
208 maybeComputeColors_();
209 }
210
212 std::size_t numDofs() const
213 { return gridGeometry_->numDofs(); }
214
216 const Problem& problem() const
217 { return *problem_; }
218
221 { return *gridGeometry_; }
222
224 const GridView& gridView() const
225 { return gridGeometry().gridView(); }
226
229 { return *gridVariables_; }
230
233 { return *gridVariables_; }
234
237 { return *jacobian_; }
238
241 { return *residual_; }
242
244 const SolutionVector& prevSol() const
245 { return *prevSol_; }
246
251 { return { LocalResidual{problem_.get(), nullptr} }; }
252
257 { gridVariables().update(cursol); }
258
262 void resetTimeStep(const SolutionVector &cursol)
263 {
264 gridVariables().resetTimeStep(cursol);
265 this->clearStages();
266 }
267
269 {
270 spatialOperatorEvaluations_.clear();
271 temporalOperatorEvaluations_.clear();
272 stageParams_.reset();
273 }
274
275 template<class StageParams>
277 {
278 stageParams_ = std::move(params);
279 const auto curStage = stageParams_->size() - 1;
280
281 // in the first stage, also assemble the residual
282 // at the previous time level (stage 0 residual)
283 if (curStage == 1)
284 {
285 // update time in variables?
286 setProblemTime_(*problem_, stageParams_->timeAtStage(0));
287
288 resetResidual_(); // residual resized and zero
289
290 assert(spatialOperatorEvaluations_.size() >= 0);
291 if (spatialOperatorEvaluations_.size() == 0)
292 {
293 spatialOperatorEvaluations_.push_back(*residual_);
294 temporalOperatorEvaluations_.push_back(*residual_);
295
296 // assemble stage 0 residuals
297 assemble_([&](const auto& element)
298 {
299 LocalAssembler localAssembler(*this, element, *prevSol_);
300 localAssembler.localResidual().spatialWeight(1.0);
301 localAssembler.localResidual().temporalWeight(1.0);
302 localAssembler.assembleCurrentResidual(spatialOperatorEvaluations_.back(), temporalOperatorEvaluations_.back());
303 });
304 }
305
306 // we don't delete the first stage so it can be reused in a restarted
307 // time integration step. The evaluations are only deleted
308 // when explicitly requested by calling clearStages().
309 // So if here the vector is non-empty, we don't need to evaluate again
310 // (this should only occur if we are restarting time integration, e.g.
311 // with a different time step size)
312 else if (spatialOperatorEvaluations_.size() > 0)
313 {
315 spatialOperatorEvaluations_.resize(1);
316 temporalOperatorEvaluations_.resize(1);
317 }
318 }
319
320 // update time in variables?
321 setProblemTime_(*problem_, stageParams_->timeAtStage(curStage));
322
323 resetResidual_(); // residual resized and zero
324
325 if (spatialOperatorEvaluations_.size() != curStage)
326 DUNE_THROW(Dune::InvalidStateException,
327 "Invalid state. Maybe you forgot to call clearStages()");
328
329 // allocate memory for this stage
330 spatialOperatorEvaluations_.push_back(*residual_);
331 temporalOperatorEvaluations_.push_back(*residual_);
332 }
333
336 { return false; }
337
338 bool isImplicit() const
339 { return timeSteppingMethod_->implicit(); }
340
341private:
345 void setJacobianPattern_()
346 {
347 // resize the jacobian and the residual
348 const auto numDofs = this->numDofs();
349 jacobian_->setSize(numDofs, numDofs);
350
351 // create occupation pattern of the jacobian
352 if (timeSteppingMethod_->implicit())
353 getJacobianPattern<true>(gridGeometry()).exportIdx(*jacobian_);
354 else
355 getJacobianPattern<false>(gridGeometry()).exportIdx(*jacobian_);
356 }
357
359 void setResidualSize_(ResidualType& res)
360 { res.resize(numDofs()); }
361
363 void maybeComputeColors_()
364 {
365 if (enableMultithreading_)
366 elementSets_ = computeColoring(gridGeometry()).sets;
367 }
368
369 // reset the residual vector to 0.0
370 void resetResidual_()
371 {
372 if(!residual_)
373 {
374 residual_ = std::make_shared<ResidualType>();
375 setResidualSize_(*residual_);
376 }
377
378 setResidualSize_(constrainedDofs_);
379
380 (*residual_) = 0.0;
381 constrainedDofs_ = 0.0;
382 }
383
384 // reset the Jacobian matrix to 0.0
385 template <class PartialReassembler = DefaultPartialReassembler>
386 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
387 {
388 if(!jacobian_)
389 {
390 jacobian_ = std::make_shared<JacobianMatrix>();
391 jacobian_->setBuildMode(JacobianMatrix::random);
392 setJacobianPattern_();
393 }
394
395 if (partialReassembler)
396 partialReassembler->resetJacobian(*this);
397 else
398 *jacobian_ = 0.0;
399 }
400
406 template<typename AssembleElementFunc>
407 void assemble_(AssembleElementFunc&& assembleElement) const
408 {
409 // a state that will be checked on all processes
410 bool succeeded = false;
411
412 // try assembling using the local assembly function
413 try
414 {
415 if (enableMultithreading_)
416 {
417 assert(elementSets_.size() > 0);
418
419 // make this element loop run in parallel
420 // for this we have to color the elements so that we don't get
421 // race conditions when writing into the global matrix
422 // each color can be assembled using multiple threads
423 for (const auto& elements : elementSets_)
424 {
425 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
426 {
427 const auto element = gridView().grid().entity(elements[i]);
428 assembleElement(element);
429 });
430 }
431 }
432 else
433 for (const auto& element : elements(gridView()))
434 assembleElement(element);
435
436 // if we get here, everything worked well on this process
437 succeeded = true;
438 }
439 // throw exception if a problem occurred
440 catch (NumericalProblem &e)
441 {
442 std::cout << "rank " << gridView().comm().rank()
443 << " caught an exception while assembling:" << e.what()
444 << "\n";
445 succeeded = false;
446 }
447
448 // make sure everything worked well on all processes
449 if (gridView().comm().size() > 1)
450 succeeded = gridView().comm().min(succeeded);
451
452 // if not succeeded rethrow the error on all processes
453 if (!succeeded)
454 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
455 }
456
457 // TODO make this nicer with a is_detected trait in a common location
458 template<class P>
459 void setProblemTime_(const P& p, const Scalar t)
460 { setProblemTimeImpl_(p, t, 0); }
461
462 template<class P>
463 auto setProblemTimeImpl_(const P& p, const Scalar t, int) -> decltype(p.setTime(0))
464 { p.setTime(t); }
465
466 template<class P>
467 void setProblemTimeImpl_(const P& p, const Scalar t, long)
468 {}
469
470 std::shared_ptr<const Experimental::MultiStageMethod<Scalar>> timeSteppingMethod_;
471 std::vector<ResidualType> spatialOperatorEvaluations_;
472 std::vector<ResidualType> temporalOperatorEvaluations_;
473 ResidualType constrainedDofs_;
474 std::shared_ptr<const StageParams> stageParams_;
475
477 std::shared_ptr<const Problem> problem_;
478
480 std::shared_ptr<const GridGeometry> gridGeometry_;
481
483 std::shared_ptr<GridVariables> gridVariables_;
484
486 const SolutionVector* prevSol_ = nullptr;
487
489 std::shared_ptr<JacobianMatrix> jacobian_;
490 std::shared_ptr<ResidualType> residual_;
491
493 bool enableMultithreading_ = false;
494 std::deque<std::vector<ElementSeed>> elementSets_;
495};
496
497} // namespace Dumux
498
499#endif
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: experimental/assembly/cclocalassembler.hh:181
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: experimental/assembly/cvfelocalassembler.hh:326
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multistagefvassembler.hh:89
const Problem & problem() const
The problem.
Definition: multistagefvassembler.hh:216
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multistagefvassembler.hh:184
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: multistagefvassembler.hh:212
MultiStageFVLocalOperator< LocalResidual > localResidual() const
Create a local residual object (used by the local assembler)
Definition: multistagefvassembler.hh:250
typename Dumux::Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: multistagefvassembler.hh:106
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: multistagefvassembler.hh:262
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: multistagefvassembler.hh:104
ResidualType & residual()
The residual vector (rhs)
Definition: multistagefvassembler.hh:240
void clearStages()
Definition: multistagefvassembler.hh:268
MultiStageFVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const Experimental::MultiStageMethod< Scalar > > msMethod, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: multistagefvassembler.hh:118
bool isStationaryProblem() const
TODO get rid of this (called by Newton but shouldn't be necessary)
Definition: multistagefvassembler.hh:335
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multistagefvassembler.hh:191
void prepareStage(SolutionVector &x, StageParams params)
Definition: multistagefvassembler.hh:276
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: multistagefvassembler.hh:108
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multistagefvassembler.hh:204
const GridView & gridView() const
The gridview.
Definition: multistagefvassembler.hh:224
GridGeo GridGeometry
Definition: multistagefvassembler.hh:110
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: multistagefvassembler.hh:244
GridVariables & gridVariables()
The global grid variables.
Definition: multistagefvassembler.hh:228
GetPropType< TypeTag, Properties::Problem > Problem
Definition: multistagefvassembler.hh:111
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: multistagefvassembler.hh:256
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: multistagefvassembler.hh:102
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: multistagefvassembler.hh:105
bool isImplicit() const
Definition: multistagefvassembler.hh:338
const GridVariables & gridVariables() const
The global grid variables.
Definition: multistagefvassembler.hh:232
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multistagefvassembler.hh:142
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: multistagefvassembler.hh:220
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: multistagefvassembler.hh:236
Definition: multistagefvlocaloperator.hh:23
Abstract interface for one-step multi-stage method parameters in Shu/Osher form.
Definition: multistagemethods.hh:75
Data object for the parameters of a given stage.
Definition: multistagetimestepper.hh:31
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper to extract native Dune vector types from particular Dumux types.
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An assembler for Jacobian and residual contribution per element (CVFE methods)
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
A local operator wrapper for multi-stage time stepping schemes.
Parameters for different multistage time stepping methods.
A time stepper performing a single time step of a transient simulation.
Multithreading in Dumux.
constexpr Box box
Definition: method.hh:147
Definition: multistagefvassembler.hh:45
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod > LocalAssemblerChooser_t
Definition: multistagefvassembler.hh:74
Definition: experimental/assembly/cclocalassembler.hh:36
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:74
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:239
Parallel for loop (multithreading)
Provides a helper class for nonoverlapping decomposition.
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Definition: multistagefvassembler.hh:48
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:292
Type traits to be used with vector types.