version 3.9
assembly/fvassembler.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//
12#ifndef DUMUX_FV_ASSEMBLER_HH
13#define DUMUX_FV_ASSEMBLER_HH
14
15#include <vector>
16#include <deque>
17#include <type_traits>
18#include <memory>
19
20#include <dune/istl/matrixindexset.hh>
21
27
31
35
38
39#include "cvfelocalassembler.hh"
40#include "cclocalassembler.hh"
41#include "fclocalassembler.hh"
42
43namespace Dumux::Detail {
44
45template<class DiscretizationMethod>
47
48template<class DM>
49struct LocalAssemblerChooser<DiscretizationMethods::CVFE<DM>>
50{
51 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
53};
54
55template<>
56struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa>
57{
58 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
60};
61
62template<>
63struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa>
64{
65 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
67};
68
69template<>
70struct LocalAssemblerChooser<DiscretizationMethods::FCStaggered>
71{
72 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
74};
75
76template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
79>::template type<TypeTag, Impl, diffMethod, isImplicit>;
80
81} // end namespace Dumux::Detail
82
83namespace Dumux {
84
92template<class TypeTag, DiffMethod diffMethod, bool isImplicit = true>
94{
96 using GridView = typename GridGeo::GridView;
98 using Element = typename GridView::template Codim<0>::Entity;
99 using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed;
101
102 static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box;
103
106
107public:
112
114
115 using GridGeometry = GridGeo;
117
123 FVAssembler(std::shared_ptr<const Problem> problem,
124 std::shared_ptr<const GridGeometry> gridGeometry,
125 std::shared_ptr<GridVariables> gridVariables)
126 : problem_(problem)
127 , gridGeometry_(gridGeometry)
128 , gridVariables_(gridVariables)
129 , timeLoop_()
130 , isStationaryProblem_(true)
131 {
132 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
134 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
136 && getParam<bool>("Assembly.Multithreading", true);
137
138 maybeComputeColors_();
139 }
140
146 FVAssembler(std::shared_ptr<const Problem> problem,
147 std::shared_ptr<const GridGeometry> gridGeometry,
148 std::shared_ptr<GridVariables> gridVariables,
149 std::shared_ptr<const TimeLoop> timeLoop,
150 const SolutionVector& prevSol)
151 : problem_(problem)
152 , gridGeometry_(gridGeometry)
153 , gridVariables_(gridVariables)
154 , timeLoop_(timeLoop)
155 , prevSol_(&prevSol)
156 , isStationaryProblem_(!timeLoop)
157 {
159 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
161 && getParam<bool>("Assembly.Multithreading", true);
162
163 maybeComputeColors_();
164 }
165
170 template<class PartialReassembler = DefaultPartialReassembler>
171 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
172 {
173 checkAssemblerState_();
174 resetJacobian_(partialReassembler);
175 resetResidual_();
176
177 assemble_([&](const Element& element)
178 {
179 LocalAssembler localAssembler(*this, element, curSol);
180 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
181 });
182
183 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
184 }
185
190 {
191 checkAssemblerState_();
192 resetJacobian_();
193
194 assemble_([&](const Element& element)
195 {
196 LocalAssembler localAssembler(*this, element, curSol);
197 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
198 });
199 }
200
203 {
204 resetResidual_();
205 assembleResidual(*residual_, curSol);
206 }
207
209 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
210 {
211 checkAssemblerState_();
212
213 assemble_([&](const Element& element)
214 {
215 LocalAssembler localAssembler(*this, element, curSol);
216 localAssembler.assembleResidual(r);
217 });
218 }
219
225 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
226 std::shared_ptr<ResidualType> r)
227 {
228 jacobian_ = A;
229 residual_ = r;
230
231 // check and/or set the BCRS matrix's build mode
232 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
233 jacobian_->setBuildMode(JacobianMatrix::random);
234 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
235 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
236
237 setResidualSize_();
238 setJacobianPattern_();
239 }
240
246 {
247 jacobian_ = std::make_shared<JacobianMatrix>();
248 jacobian_->setBuildMode(JacobianMatrix::random);
249 residual_ = std::make_shared<ResidualType>();
250
251 setResidualSize_();
252 setJacobianPattern_();
253 }
254
259 {
260 setResidualSize_();
261 setJacobianPattern_();
262 maybeComputeColors_();
263 }
264
266 std::size_t numDofs() const
267 { return gridGeometry_->numDofs(); }
268
270 const Problem& problem() const
271 { return *problem_; }
272
275 { return *gridGeometry_; }
276
278 const GridView& gridView() const
279 { return gridGeometry().gridView(); }
280
283 { return *gridVariables_; }
284
287 { return *gridVariables_; }
288
291 { return *jacobian_; }
292
295 { return *residual_; }
296
298 const SolutionVector& prevSol() const
299 { return *prevSol_; }
300
305 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
306 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
307
313 { prevSol_ = &u; }
314
319 { return isStationaryProblem_; }
320
324 LocalResidual localResidual() const
325 { return LocalResidual(problem_.get(), timeLoop_.get()); }
326
331 {
332 gridVariables().update(cursol);
333 }
334
338 void resetTimeStep(const SolutionVector &cursol)
339 {
340 gridVariables().resetTimeStep(cursol);
341 }
342
343private:
347 void setJacobianPattern_()
348 {
349 // resize the jacobian and the residual
350 const auto numDofs = this->numDofs();
351 jacobian_->setSize(numDofs, numDofs);
352
353 // create occupation pattern of the jacobian
354 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
355
356 // export pattern to jacobian
357 occupationPattern.exportIdx(*jacobian_);
358 }
359
361 void setResidualSize_()
362 { residual_->resize(numDofs()); }
363
365 void maybeComputeColors_()
366 {
367 if (enableMultithreading_)
368 elementSets_ = computeColoring(gridGeometry()).sets;
369 }
370
371 // reset the residual vector to 0.0
372 void resetResidual_()
373 {
374 if(!residual_)
375 {
376 residual_ = std::make_shared<ResidualType>();
377 setResidualSize_();
378 }
379
380 (*residual_) = 0.0;
381 }
382
383 // reset the Jacobian matrix to 0.0
384 template <class PartialReassembler = DefaultPartialReassembler>
385 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
386 {
387 if(!jacobian_)
388 {
389 jacobian_ = std::make_shared<JacobianMatrix>();
390 jacobian_->setBuildMode(JacobianMatrix::random);
391 setJacobianPattern_();
392 }
393
394 if (partialReassembler)
395 partialReassembler->resetJacobian(*this);
396 else
397 *jacobian_ = 0.0;
398 }
399
400 // check if the assembler is in a correct state for assembly
401 void checkAssemblerState_() const
402 {
403 if (!isStationaryProblem_ && !prevSol_)
404 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
405 }
406
412 template<typename AssembleElementFunc>
413 void assemble_(AssembleElementFunc&& assembleElement) const
414 {
415 // a state that will be checked on all processes
416 bool succeeded = false;
417
418 // try assembling using the local assembly function
419 try
420 {
421 if (enableMultithreading_)
422 {
423 assert(elementSets_.size() > 0);
424
425 // make this element loop run in parallel
426 // for this we have to color the elements so that we don't get
427 // race conditions when writing into the global matrix
428 // each color can be assembled using multiple threads
429 for (const auto& elements : elementSets_)
430 {
431 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
432 {
433 const auto element = gridView().grid().entity(elements[i]);
434 assembleElement(element);
435 });
436 }
437 }
438 else
439 for (const auto& element : elements(gridView()))
440 assembleElement(element);
441
442 // if we get here, everything worked well on this process
443 succeeded = true;
444 }
445 // throw exception if a problem occurred
446 catch (NumericalProblem &e)
447 {
448 std::cout << "rank " << gridView().comm().rank()
449 << " caught an exception while assembling:" << e.what()
450 << "\n";
451 succeeded = false;
452 }
453
454 // make sure everything worked well on all processes
455 if (gridView().comm().size() > 1)
456 succeeded = gridView().comm().min(succeeded);
457
458 // if not succeeded rethrow the error on all processes
459 if (!succeeded)
460 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
461 }
462
463 template<class GG>
464 void enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry)
465 {
466 if constexpr (Detail::hasPeriodicDofMap<GG>())
467 {
468 for (const auto& m : gridGeometry.periodicDofMap())
469 {
470 if (m.first < m.second)
471 {
472 // add the second row to the first
473 res[m.first] += res[m.second];
474 const auto end = jac[m.second].end();
475 for (auto it = jac[m.second].begin(); it != end; ++it)
476 jac[m.first][it.index()] += (*it);
477
478 // enforce constraint in second row
479 res[m.second] = curSol[m.second] - curSol[m.first];
480
481 // set derivatives accordingly in jacobian, i.e. id for m.second and -id for m.first
482 auto setMatrixBlock = [] (auto& matrixBlock, double diagValue)
483 {
484 for (int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
485 matrixBlock[eIdx][eIdx] = diagValue;
486 };
487
488 for (auto it = jac[m.second].begin(); it != end; ++it)
489 {
490 auto& matrixBlock = *it;
491 matrixBlock = 0.0;
492
493 assert(matrixBlock.N() == matrixBlock.M());
494 if(it.index() == m.second)
495 setMatrixBlock(matrixBlock, 1.0);
496
497 if(it.index() == m.first)
498 setMatrixBlock(matrixBlock, -1.0);
499
500 }
501 }
502 }
503 }
504 }
505
507 std::shared_ptr<const Problem> problem_;
508
510 std::shared_ptr<const GridGeometry> gridGeometry_;
511
513 std::shared_ptr<GridVariables> gridVariables_;
514
516 std::shared_ptr<const TimeLoop> timeLoop_;
517
519 const SolutionVector* prevSol_ = nullptr;
520
522 bool isStationaryProblem_;
523
525 std::shared_ptr<JacobianMatrix> jacobian_;
526 std::shared_ptr<ResidualType> residual_;
527
529 bool enableMultithreading_ = false;
530 std::deque<std::vector<ElementSeed>> elementSets_;
531};
532
533} // namespace Dumux
534
535#endif
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: assembly/cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:302
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:94
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:294
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:108
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:189
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/fvassembler.hh:258
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:278
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:270
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: assembly/fvassembler.hh:111
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:245
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:318
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:109
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:324
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/fvassembler.hh:110
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: assembly/fvassembler.hh:146
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:266
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< ResidualType > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: assembly/fvassembler.hh:225
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:282
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:290
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:209
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:116
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:305
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:202
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:113
void setPreviousSolution(const SolutionVector &u)
Sets the solution from which to start the time integration. Has to be called prior to assembly for ti...
Definition: assembly/fvassembler.hh:312
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:274
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:338
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:286
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: assembly/fvassembler.hh:171
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:115
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables)
The constructor for stationary problems.
Definition: assembly/fvassembler.hh:123
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:298
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:330
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:284
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:84
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
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)
An assembler for Jacobian and residual contribution per element (face-centered staggered 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.
Multithreading in Dumux.
Distance implementation details.
Definition: cvfelocalresidual.hh:25
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition: assembly/fvassembler.hh:79
constexpr Box box
Definition: method.hh:147
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:74
Definition: adapt.hh:17
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.
Type traits to detect periodicity support.
Definition: assembly/fvassembler.hh:46
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:292
Type traits to be used with vector types.