version 3.8
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 { gridVariables().resetTimeStep(cursol); }
264
266 {
267 spatialOperatorEvaluations_.clear();
268 temporalOperatorEvaluations_.clear();
269 stageParams_.reset();
270 }
271
272 template<class StageParams>
274 {
275 stageParams_ = std::move(params);
276 const auto curStage = stageParams_->size() - 1;
277
278 // in the first stage, also assemble the old residual
279 if (curStage == 1)
280 {
281 // update time in variables?
282 setProblemTime_(*problem_, stageParams_->timeAtStage(curStage));
283
284 resetResidual_(); // residual resized and zero
285 spatialOperatorEvaluations_.push_back(*residual_);
286 temporalOperatorEvaluations_.push_back(*residual_);
287
288 // assemble stage 0 residuals
289 assemble_([&](const auto& element)
290 {
291 LocalAssembler localAssembler(*this, element, *prevSol_);
292 localAssembler.localResidual().spatialWeight(1.0);
293 localAssembler.localResidual().temporalWeight(1.0);
294 localAssembler.assembleCurrentResidual(spatialOperatorEvaluations_.back(), temporalOperatorEvaluations_.back());
295 });
296 }
297
298 // update time in variables?
299 setProblemTime_(*problem_, stageParams_->timeAtStage(curStage));
300
301 resetResidual_(); // residual resized and zero
302 spatialOperatorEvaluations_.push_back(*residual_);
303 temporalOperatorEvaluations_.push_back(*residual_);
304 }
305
308 { return false; }
309
310 bool isImplicit() const
311 { return timeSteppingMethod_->implicit(); }
312
313private:
317 void setJacobianPattern_()
318 {
319 // resize the jacobian and the residual
320 const auto numDofs = this->numDofs();
321 jacobian_->setSize(numDofs, numDofs);
322
323 // create occupation pattern of the jacobian
324 if (timeSteppingMethod_->implicit())
325 getJacobianPattern<true>(gridGeometry()).exportIdx(*jacobian_);
326 else
327 getJacobianPattern<false>(gridGeometry()).exportIdx(*jacobian_);
328 }
329
331 void setResidualSize_(ResidualType& res)
332 { res.resize(numDofs()); }
333
335 void maybeComputeColors_()
336 {
337 if (enableMultithreading_)
338 elementSets_ = computeColoring(gridGeometry()).sets;
339 }
340
341 // reset the residual vector to 0.0
342 void resetResidual_()
343 {
344 if(!residual_)
345 {
346 residual_ = std::make_shared<ResidualType>();
347 setResidualSize_(*residual_);
348 }
349
350 setResidualSize_(constrainedDofs_);
351
352 (*residual_) = 0.0;
353 constrainedDofs_ = 0.0;
354 }
355
356 // reset the Jacobian matrix to 0.0
357 template <class PartialReassembler = DefaultPartialReassembler>
358 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
359 {
360 if(!jacobian_)
361 {
362 jacobian_ = std::make_shared<JacobianMatrix>();
363 jacobian_->setBuildMode(JacobianMatrix::random);
364 setJacobianPattern_();
365 }
366
367 if (partialReassembler)
368 partialReassembler->resetJacobian(*this);
369 else
370 *jacobian_ = 0.0;
371 }
372
378 template<typename AssembleElementFunc>
379 void assemble_(AssembleElementFunc&& assembleElement) const
380 {
381 // a state that will be checked on all processes
382 bool succeeded = false;
383
384 // try assembling using the local assembly function
385 try
386 {
387 if (enableMultithreading_)
388 {
389 assert(elementSets_.size() > 0);
390
391 // make this element loop run in parallel
392 // for this we have to color the elements so that we don't get
393 // race conditions when writing into the global matrix
394 // each color can be assembled using multiple threads
395 for (const auto& elements : elementSets_)
396 {
397 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
398 {
399 const auto element = gridView().grid().entity(elements[i]);
400 assembleElement(element);
401 });
402 }
403 }
404 else
405 for (const auto& element : elements(gridView()))
406 assembleElement(element);
407
408 // if we get here, everything worked well on this process
409 succeeded = true;
410 }
411 // throw exception if a problem occurred
412 catch (NumericalProblem &e)
413 {
414 std::cout << "rank " << gridView().comm().rank()
415 << " caught an exception while assembling:" << e.what()
416 << "\n";
417 succeeded = false;
418 }
419
420 // make sure everything worked well on all processes
421 if (gridView().comm().size() > 1)
422 succeeded = gridView().comm().min(succeeded);
423
424 // if not succeeded rethrow the error on all processes
425 if (!succeeded)
426 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
427 }
428
429 // TODO make this nicer with a is_detected trait in a common location
430 template<class P>
431 void setProblemTime_(const P& p, const Scalar t)
432 { setProblemTimeImpl_(p, t, 0); }
433
434 template<class P>
435 auto setProblemTimeImpl_(const P& p, const Scalar t, int) -> decltype(p.setTime(0))
436 { p.setTime(t); }
437
438 template<class P>
439 void setProblemTimeImpl_(const P& p, const Scalar t, long)
440 {}
441
442 std::shared_ptr<const Experimental::MultiStageMethod<Scalar>> timeSteppingMethod_;
443 std::vector<ResidualType> spatialOperatorEvaluations_;
444 std::vector<ResidualType> temporalOperatorEvaluations_;
445 ResidualType constrainedDofs_;
446 std::shared_ptr<const StageParams> stageParams_;
447
449 std::shared_ptr<const Problem> problem_;
450
452 std::shared_ptr<const GridGeometry> gridGeometry_;
453
455 std::shared_ptr<GridVariables> gridVariables_;
456
458 const SolutionVector* prevSol_ = nullptr;
459
461 std::shared_ptr<JacobianMatrix> jacobian_;
462 std::shared_ptr<ResidualType> residual_;
463
465 bool enableMultithreading_ = false;
466 std::deque<std::vector<ElementSeed>> elementSets_;
467};
468
469} // namespace Dumux
470
471#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:265
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:307
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:273
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:310
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:29
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.