3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FV_ASSEMBLER_HH
25#define DUMUX_FV_ASSEMBLER_HH
26
27#include <vector>
28#include <deque>
29#include <type_traits>
30#include <memory>
31
32#include <dune/istl/matrixindexset.hh>
33
37
40
44
47
48#include "boxlocalassembler.hh"
49#include "cclocalassembler.hh"
50#include "fclocalassembler.hh"
53
54namespace Dumux::Detail {
55
56template<class DiscretizationMethod>
58
59template<>
60struct LocalAssemblerChooser<DiscretizationMethods::Box>
61{
62 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
64};
65
66template<>
67struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa>
68{
69 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
71};
72
73template<>
74struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa>
75{
76 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
78};
79
80template<>
81struct LocalAssemblerChooser<DiscretizationMethods::FCStaggered>
82{
83 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
85};
86
87template<>
88struct LocalAssemblerChooser<DiscretizationMethods::FCDiamond>
89{
90 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
92};
93
94template<>
95struct LocalAssemblerChooser<DiscretizationMethods::PQ1Bubble>
96{
97 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
99};
100
101template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
104>::template type<TypeTag, Impl, diffMethod, isImplicit>;
105
106} // end namespace Dumux::Detail
107
108namespace Dumux {
109
117template<class TypeTag, DiffMethod diffMethod, bool isImplicit = true>
119{
121 using GridView = typename GridGeo::GridView;
123 using Element = typename GridView::template Codim<0>::Entity;
124 using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed;
127
128 static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box;
129
132
133public:
136 using GridGeometry = GridGeo;
139
140 using ResidualType = SolutionVector;
141
147 FVAssembler(std::shared_ptr<const Problem> problem,
148 std::shared_ptr<const GridGeometry> gridGeometry,
149 std::shared_ptr<GridVariables> gridVariables)
150 : problem_(problem)
151 , gridGeometry_(gridGeometry)
152 , gridVariables_(gridVariables)
153 , timeLoop_()
154 , isStationaryProblem_(true)
155 {
156 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
158 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
160 && getParam<bool>("Assembly.Multithreading", true);
161
162 maybeComputeColors_();
163 }
164
170 FVAssembler(std::shared_ptr<const Problem> problem,
171 std::shared_ptr<const GridGeometry> gridGeometry,
172 std::shared_ptr<GridVariables> gridVariables,
173 std::shared_ptr<const TimeLoop> timeLoop,
174 const SolutionVector& prevSol)
175 : problem_(problem)
176 , gridGeometry_(gridGeometry)
177 , gridVariables_(gridVariables)
178 , timeLoop_(timeLoop)
179 , prevSol_(&prevSol)
180 , isStationaryProblem_(!timeLoop)
181 {
183 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
185 && getParam<bool>("Assembly.Multithreading", true);
186
187 maybeComputeColors_();
188 }
189
194 template<class PartialReassembler = DefaultPartialReassembler>
195 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
196 {
197 checkAssemblerState_();
198 resetJacobian_(partialReassembler);
199 resetResidual_();
200
201 assemble_([&](const Element& element)
202 {
203 LocalAssembler localAssembler(*this, element, curSol);
204 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
205 });
206
207 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
208 }
209
213 void assembleJacobian(const SolutionVector& curSol)
214 {
215 checkAssemblerState_();
216 resetJacobian_();
217
218 assemble_([&](const Element& element)
219 {
220 LocalAssembler localAssembler(*this, element, curSol);
221 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
222 });
223 }
224
226 void assembleResidual(const SolutionVector& curSol)
227 {
228 resetResidual_();
229 assembleResidual(*residual_, curSol);
230 }
231
233 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
234 {
235 checkAssemblerState_();
236
237 assemble_([&](const Element& element)
238 {
239 LocalAssembler localAssembler(*this, element, curSol);
240 localAssembler.assembleResidual(r);
241 });
242 }
243
245 Scalar residualNorm(const SolutionVector& curSol) const
246 {
248 assembleResidual(residual, curSol);
249
250 // issue a warning if the calculation is used in parallel with overlap
251 static bool warningIssued = false;
252
253 if (gridView().comm().size() > 1 && gridView().overlapSize(0) == 0)
254 {
255 if constexpr (isBox)
256 {
257 using DM = typename GridGeometry::VertexMapper;
259
260 PVHelper vectorHelper(gridView(), gridGeometry_->vertexMapper());
261
262 vectorHelper.makeNonOverlappingConsistent(residual);
263 }
264 }
265 else if (!warningIssued)
266 {
267 if (gridView().comm().size() > 1 && gridView().comm().rank() == 0)
268 std::cout << "\nWarning: norm calculation adds entries corresponding to\n"
269 << "overlapping entities multiple times. Please use the norm\n"
270 << "function provided by a linear solver instead." << std::endl;
271
272 warningIssued = true;
273 }
274
275 // calculate the square norm of the residual
276 Scalar result2 = residual.two_norm2();
277 if (gridView().comm().size() > 1)
278 result2 = gridView().comm().sum(result2);
279
280 using std::sqrt;
281 return sqrt(result2);
282 }
283
289 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
290 std::shared_ptr<SolutionVector> r)
291 {
292 jacobian_ = A;
293 residual_ = r;
294
295 // check and/or set the BCRS matrix's build mode
296 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
297 jacobian_->setBuildMode(JacobianMatrix::random);
298 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
299 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
300
301 setResidualSize_();
302 setJacobianPattern_();
303 }
304
310 {
311 jacobian_ = std::make_shared<JacobianMatrix>();
312 jacobian_->setBuildMode(JacobianMatrix::random);
313 residual_ = std::make_shared<SolutionVector>();
314
315 setResidualSize_();
316 setJacobianPattern_();
317 }
318
323 {
324 setResidualSize_();
325 setJacobianPattern_();
326 maybeComputeColors_();
327 }
328
330 std::size_t numDofs() const
331 { return gridGeometry_->numDofs(); }
332
334 const Problem& problem() const
335 { return *problem_; }
336
339 { return *gridGeometry_; }
340
342 const GridView& gridView() const
343 { return gridGeometry().gridView(); }
344
347 { return *gridVariables_; }
348
351 { return *gridVariables_; }
352
355 { return *jacobian_; }
356
358 SolutionVector& residual()
359 { return *residual_; }
360
362 const SolutionVector& prevSol() const
363 { return *prevSol_; }
364
369 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
370 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
371
376 void setPreviousSolution(const SolutionVector& u)
377 { prevSol_ = &u; }
378
383 { return isStationaryProblem_; }
384
388 LocalResidual localResidual() const
389 { return LocalResidual(problem_.get(), timeLoop_.get()); }
390
394 void updateGridVariables(const SolutionVector &cursol)
395 {
396 gridVariables().update(cursol);
397 }
398
402 void resetTimeStep(const SolutionVector &cursol)
403 {
404 gridVariables().resetTimeStep(cursol);
405 }
406
407private:
411 void setJacobianPattern_()
412 {
413 // resize the jacobian and the residual
414 const auto numDofs = this->numDofs();
415 jacobian_->setSize(numDofs, numDofs);
416
417 // create occupation pattern of the jacobian
418 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
419
420 // export pattern to jacobian
421 occupationPattern.exportIdx(*jacobian_);
422 }
423
425 void setResidualSize_()
426 { residual_->resize(numDofs()); }
427
429 void maybeComputeColors_()
430 {
431 if (enableMultithreading_)
432 elementSets_ = computeColoring(gridGeometry()).sets;
433 }
434
435 // reset the residual vector to 0.0
436 void resetResidual_()
437 {
438 if(!residual_)
439 {
440 residual_ = std::make_shared<SolutionVector>();
441 setResidualSize_();
442 }
443
444 (*residual_) = 0.0;
445 }
446
447 // reset the Jacobian matrix to 0.0
448 template <class PartialReassembler = DefaultPartialReassembler>
449 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
450 {
451 if(!jacobian_)
452 {
453 jacobian_ = std::make_shared<JacobianMatrix>();
454 jacobian_->setBuildMode(JacobianMatrix::random);
455 setJacobianPattern_();
456 }
457
458 if (partialReassembler)
459 partialReassembler->resetJacobian(*this);
460 else
461 *jacobian_ = 0.0;
462 }
463
464 // check if the assembler is in a correct state for assembly
465 void checkAssemblerState_() const
466 {
467 if (!isStationaryProblem_ && !prevSol_)
468 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
469 }
470
476 template<typename AssembleElementFunc>
477 void assemble_(AssembleElementFunc&& assembleElement) const
478 {
479 // a state that will be checked on all processes
480 bool succeeded = false;
481
482 // try assembling using the local assembly function
483 try
484 {
485 if (enableMultithreading_)
486 {
487 assert(elementSets_.size() > 0);
488
489 // make this element loop run in parallel
490 // for this we have to color the elements so that we don't get
491 // race conditions when writing into the global matrix
492 // each color can be assembled using multiple threads
493 for (const auto& elements : elementSets_)
494 {
495 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
496 {
497 const auto element = gridView().grid().entity(elements[i]);
498 assembleElement(element);
499 });
500 }
501 }
502 else
503 for (const auto& element : elements(gridView()))
504 assembleElement(element);
505
506 // if we get here, everything worked well on this process
507 succeeded = true;
508 }
509 // throw exception if a problem occurred
510 catch (NumericalProblem &e)
511 {
512 std::cout << "rank " << gridView().comm().rank()
513 << " caught an exception while assembling:" << e.what()
514 << "\n";
515 succeeded = false;
516 }
517
518 // make sure everything worked well on all processes
519 if (gridView().comm().size() > 1)
520 succeeded = gridView().comm().min(succeeded);
521
522 // if not succeeded rethrow the error on all processes
523 if (!succeeded)
524 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
525 }
526
527 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethods::box, void>
528 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry)
529 {
530 for (const auto& m : gridGeometry.periodicVertexMap())
531 {
532 if (m.first < m.second)
533 {
534 // add the second row to the first
535 res[m.first] += res[m.second];
536 const auto end = jac[m.second].end();
537 for (auto it = jac[m.second].begin(); it != end; ++it)
538 jac[m.first][it.index()] += (*it);
539
540 // enforce constraint in second row
541 res[m.second] = curSol[m.second] - curSol[m.first];
542 for (auto it = jac[m.second].begin(); it != end; ++it)
543 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
544 }
545 }
546 }
547
548 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethods::box, void>
549 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry) {}
550
552 std::shared_ptr<const Problem> problem_;
553
555 std::shared_ptr<const GridGeometry> gridGeometry_;
556
558 std::shared_ptr<GridVariables> gridVariables_;
559
561 std::shared_ptr<const TimeLoop> timeLoop_;
562
564 const SolutionVector* prevSol_ = nullptr;
565
567 bool isStationaryProblem_;
568
570 std::shared_ptr<JacobianMatrix> jacobian_;
571 std::shared_ptr<SolutionVector> residual_;
572
574 bool enableMultithreading_ = false;
575 std::deque<std::vector<ElementSeed>> elementSets_;
576};
577
578} // namespace Dumux
579
580#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Coloring schemes for shared-memory-parallel assembly.
An assembler for Jacobian and residual contribution per element (face-centered diamond methods)
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
Helper function to generate Jacobian pattern for different discretization methods.
An assembler for Jacobian and residual contribution per element (pq1bubble method)
dune-grid capabilities compatibility layer
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
Multithreading in Dumux.
Parallel for loop (multithreading)
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:57
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:172
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:251
Distance implementation details.
Definition: cvfelocalresidual.hh:37
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition: assembly/fvassembler.hh:104
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:86
CVFE< CVFEMethods::CR_RT > FCDiamond
Definition: method.hh:90
constexpr Box box
Definition: method.hh:136
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
CVFE< CVFEMethods::PQ1Bubble > PQ1Bubble
Definition: method.hh:97
An assembler for Jacobian and residual contribution per element (box methods)
Definition: boxlocalassembler.hh:257
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:136
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:304
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fcdiamondlocalassembler.hh:290
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:296
Definition: assembly/fvassembler.hh:57
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:119
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:134
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:213
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/fvassembler.hh:322
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:342
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:334
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:309
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:140
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:382
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:135
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:388
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:170
Scalar residualNorm(const SolutionVector &curSol) const
compute the residual and return it's vector norm
Definition: assembly/fvassembler.hh:245
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:330
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:346
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:354
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:233
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:137
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:369
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:226
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:138
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:376
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:338
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:402
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:350
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< SolutionVector > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: assembly/fvassembler.hh:289
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:195
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:136
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:147
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:358
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:362
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:394
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
An assembler for Jacobian and residual contribution per element (PQ1Bubble methods)
Definition: pq1bubblelocalassembler.hh:304
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:68
Definition: parallelhelpers.hh:473
Declares all properties used in Dumux.
Manages the handling of time dependent problems.