version 3.7
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
26
30
34
37
38#include "cvfelocalassembler.hh"
39#include "cclocalassembler.hh"
40#include "fclocalassembler.hh"
41
42namespace Dumux::Detail {
43
44template<class DiscretizationMethod>
46
47template<class DM>
48struct LocalAssemblerChooser<DiscretizationMethods::CVFE<DM>>
49{
50 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
52};
53
54template<>
55struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa>
56{
57 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
59};
60
61template<>
62struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa>
63{
64 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
66};
67
68template<>
69struct LocalAssemblerChooser<DiscretizationMethods::FCStaggered>
70{
71 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
73};
74
75template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
78>::template type<TypeTag, Impl, diffMethod, isImplicit>;
79
80} // end namespace Dumux::Detail
81
82namespace Dumux {
83
91template<class TypeTag, DiffMethod diffMethod, bool isImplicit = true>
93{
95 using GridView = typename GridGeo::GridView;
97 using Element = typename GridView::template Codim<0>::Entity;
98 using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed;
100
101 static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box;
102
105
106public:
111
113
114 using GridGeometry = GridGeo;
116
122 FVAssembler(std::shared_ptr<const Problem> problem,
123 std::shared_ptr<const GridGeometry> gridGeometry,
124 std::shared_ptr<GridVariables> gridVariables)
125 : problem_(problem)
126 , gridGeometry_(gridGeometry)
127 , gridVariables_(gridVariables)
128 , timeLoop_()
129 , isStationaryProblem_(true)
130 {
131 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
133 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
135 && getParam<bool>("Assembly.Multithreading", true);
136
137 maybeComputeColors_();
138 }
139
145 FVAssembler(std::shared_ptr<const Problem> problem,
146 std::shared_ptr<const GridGeometry> gridGeometry,
147 std::shared_ptr<GridVariables> gridVariables,
148 std::shared_ptr<const TimeLoop> timeLoop,
149 const SolutionVector& prevSol)
150 : problem_(problem)
151 , gridGeometry_(gridGeometry)
152 , gridVariables_(gridVariables)
153 , timeLoop_(timeLoop)
154 , prevSol_(&prevSol)
155 , isStationaryProblem_(!timeLoop)
156 {
158 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
160 && getParam<bool>("Assembly.Multithreading", true);
161
162 maybeComputeColors_();
163 }
164
169 template<class PartialReassembler = DefaultPartialReassembler>
170 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
171 {
172 checkAssemblerState_();
173 resetJacobian_(partialReassembler);
174 resetResidual_();
175
176 assemble_([&](const Element& element)
177 {
178 LocalAssembler localAssembler(*this, element, curSol);
179 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
180 });
181
182 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
183 }
184
189 {
190 checkAssemblerState_();
191 resetJacobian_();
192
193 assemble_([&](const Element& element)
194 {
195 LocalAssembler localAssembler(*this, element, curSol);
196 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
197 });
198 }
199
202 {
203 resetResidual_();
204 assembleResidual(*residual_, curSol);
205 }
206
208 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
209 {
210 checkAssemblerState_();
211
212 assemble_([&](const Element& element)
213 {
214 LocalAssembler localAssembler(*this, element, curSol);
215 localAssembler.assembleResidual(r);
216 });
217 }
218
220 [[deprecated("Use the linear solver's norm. Will be deleted after 3.7")]]
222 {
223 // issue a warning if the calculation is used in parallel with overlap
224 static bool warningIssued = false;
225
226 if (gridView().comm().size() > 1 && gridView().overlapSize(0) == 0)
227 {
228 if constexpr (isBox)
229 {
230 using DM = typename GridGeometry::VertexMapper;
232
233 PVHelper vectorHelper(gridView(), gridGeometry_->vertexMapper());
234
235 vectorHelper.makeNonOverlappingConsistent(residual);
236 }
237 }
238 else if (!warningIssued)
239 {
240 if (gridView().comm().size() > 1 && gridView().comm().rank() == 0)
241 std::cout << "\nWarning: norm calculation adds entries corresponding to\n"
242 << "overlapping entities multiple times. Please use the norm\n"
243 << "function provided by a linear solver instead." << std::endl;
244
245 warningIssued = true;
246 }
247
248 // calculate the square norm of the residual
249 Scalar result2 = residual.two_norm2();
250 if (gridView().comm().size() > 1)
251 result2 = gridView().comm().sum(result2);
252
253 using std::sqrt;
254 return sqrt(result2);
255 }
256
258 [[deprecated("Use assembleResidual and the linear solver's norm. Will be deleted after 3.7")]]
259 Scalar residualNorm(const SolutionVector& curSol) const
260 {
262 assembleResidual(residual, curSol);
263 return normOfResidual(residual);
264 }
265
271 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
272 std::shared_ptr<ResidualType> r)
273 {
274 jacobian_ = A;
275 residual_ = r;
276
277 // check and/or set the BCRS matrix's build mode
278 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
279 jacobian_->setBuildMode(JacobianMatrix::random);
280 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
281 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
282
283 setResidualSize_();
284 setJacobianPattern_();
285 }
286
292 {
293 jacobian_ = std::make_shared<JacobianMatrix>();
294 jacobian_->setBuildMode(JacobianMatrix::random);
295 residual_ = std::make_shared<ResidualType>();
296
297 setResidualSize_();
298 setJacobianPattern_();
299 }
300
305 {
306 setResidualSize_();
307 setJacobianPattern_();
308 maybeComputeColors_();
309 }
310
312 std::size_t numDofs() const
313 { return gridGeometry_->numDofs(); }
314
316 const Problem& problem() const
317 { return *problem_; }
318
321 { return *gridGeometry_; }
322
324 const GridView& gridView() const
325 { return gridGeometry().gridView(); }
326
329 { return *gridVariables_; }
330
333 { return *gridVariables_; }
334
337 { return *jacobian_; }
338
341 { return *residual_; }
342
344 const SolutionVector& prevSol() const
345 { return *prevSol_; }
346
351 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
352 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
353
359 { prevSol_ = &u; }
360
365 { return isStationaryProblem_; }
366
370 LocalResidual localResidual() const
371 { return LocalResidual(problem_.get(), timeLoop_.get()); }
372
377 {
378 gridVariables().update(cursol);
379 }
380
384 void resetTimeStep(const SolutionVector &cursol)
385 {
386 gridVariables().resetTimeStep(cursol);
387 }
388
389private:
393 void setJacobianPattern_()
394 {
395 // resize the jacobian and the residual
396 const auto numDofs = this->numDofs();
397 jacobian_->setSize(numDofs, numDofs);
398
399 // create occupation pattern of the jacobian
400 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
401
402 // export pattern to jacobian
403 occupationPattern.exportIdx(*jacobian_);
404 }
405
407 void setResidualSize_()
408 { residual_->resize(numDofs()); }
409
411 void maybeComputeColors_()
412 {
413 if (enableMultithreading_)
414 elementSets_ = computeColoring(gridGeometry()).sets;
415 }
416
417 // reset the residual vector to 0.0
418 void resetResidual_()
419 {
420 if(!residual_)
421 {
422 residual_ = std::make_shared<ResidualType>();
423 setResidualSize_();
424 }
425
426 (*residual_) = 0.0;
427 }
428
429 // reset the Jacobian matrix to 0.0
430 template <class PartialReassembler = DefaultPartialReassembler>
431 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
432 {
433 if(!jacobian_)
434 {
435 jacobian_ = std::make_shared<JacobianMatrix>();
436 jacobian_->setBuildMode(JacobianMatrix::random);
437 setJacobianPattern_();
438 }
439
440 if (partialReassembler)
441 partialReassembler->resetJacobian(*this);
442 else
443 *jacobian_ = 0.0;
444 }
445
446 // check if the assembler is in a correct state for assembly
447 void checkAssemblerState_() const
448 {
449 if (!isStationaryProblem_ && !prevSol_)
450 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
451 }
452
458 template<typename AssembleElementFunc>
459 void assemble_(AssembleElementFunc&& assembleElement) const
460 {
461 // a state that will be checked on all processes
462 bool succeeded = false;
463
464 // try assembling using the local assembly function
465 try
466 {
467 if (enableMultithreading_)
468 {
469 assert(elementSets_.size() > 0);
470
471 // make this element loop run in parallel
472 // for this we have to color the elements so that we don't get
473 // race conditions when writing into the global matrix
474 // each color can be assembled using multiple threads
475 for (const auto& elements : elementSets_)
476 {
477 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
478 {
479 const auto element = gridView().grid().entity(elements[i]);
480 assembleElement(element);
481 });
482 }
483 }
484 else
485 for (const auto& element : elements(gridView()))
486 assembleElement(element);
487
488 // if we get here, everything worked well on this process
489 succeeded = true;
490 }
491 // throw exception if a problem occurred
492 catch (NumericalProblem &e)
493 {
494 std::cout << "rank " << gridView().comm().rank()
495 << " caught an exception while assembling:" << e.what()
496 << "\n";
497 succeeded = false;
498 }
499
500 // make sure everything worked well on all processes
501 if (gridView().comm().size() > 1)
502 succeeded = gridView().comm().min(succeeded);
503
504 // if not succeeded rethrow the error on all processes
505 if (!succeeded)
506 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
507 }
508
509 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethods::box, void>
510 enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry)
511 {
512 for (const auto& m : gridGeometry.periodicVertexMap())
513 {
514 if (m.first < m.second)
515 {
516 // add the second row to the first
517 res[m.first] += res[m.second];
518 const auto end = jac[m.second].end();
519 for (auto it = jac[m.second].begin(); it != end; ++it)
520 jac[m.first][it.index()] += (*it);
521
522 // enforce constraint in second row
523 res[m.second] = curSol[m.second] - curSol[m.first];
524 for (auto it = jac[m.second].begin(); it != end; ++it)
525 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
526 }
527 }
528 }
529
530 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethods::box, void>
531 enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry) {}
532
534 std::shared_ptr<const Problem> problem_;
535
537 std::shared_ptr<const GridGeometry> gridGeometry_;
538
540 std::shared_ptr<GridVariables> gridVariables_;
541
543 std::shared_ptr<const TimeLoop> timeLoop_;
544
546 const SolutionVector* prevSol_ = nullptr;
547
549 bool isStationaryProblem_;
550
552 std::shared_ptr<JacobianMatrix> jacobian_;
553 std::shared_ptr<ResidualType> residual_;
554
556 bool enableMultithreading_ = false;
557 std::deque<std::vector<ElementSeed>> elementSets_;
558};
559
560} // namespace Dumux
561
562#endif
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: cvfelocalassembler.hh:299
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:93
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:340
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:107
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:188
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/fvassembler.hh:304
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:324
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:316
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: assembly/fvassembler.hh:110
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:291
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:364
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:108
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:370
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/fvassembler.hh:109
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:145
Scalar normOfResidual(ResidualType &residual) const
compute a residual's vector norm (this is a temporary interface introduced during the deprecation per...
Definition: assembly/fvassembler.hh:221
Scalar residualNorm(const SolutionVector &curSol) const
compute the residual and return it's vector norm
Definition: assembly/fvassembler.hh:259
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:312
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:271
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:328
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:336
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:208
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:115
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:351
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:201
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:112
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:358
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:320
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:384
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:332
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:170
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:114
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:122
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:344
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:376
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:284
Definition: parallelhelpers.hh:476
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:56
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
An assembler for Jacobian and residual contribution per element (CVFE methods)
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 (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:267
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:78
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.
Definition: assembly/fvassembler.hh:45
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.