version 3.8.0
Loading...
Searching...
No Matches
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
224 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
225 std::shared_ptr<ResidualType> r)
226 {
227 jacobian_ = A;
228 residual_ = r;
229
230 // check and/or set the BCRS matrix's build mode
231 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
232 jacobian_->setBuildMode(JacobianMatrix::random);
233 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
234 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
235
236 setResidualSize_();
237 setJacobianPattern_();
238 }
239
245 {
246 jacobian_ = std::make_shared<JacobianMatrix>();
247 jacobian_->setBuildMode(JacobianMatrix::random);
248 residual_ = std::make_shared<ResidualType>();
249
250 setResidualSize_();
251 setJacobianPattern_();
252 }
253
258 {
259 setResidualSize_();
260 setJacobianPattern_();
261 maybeComputeColors_();
262 }
263
265 std::size_t numDofs() const
266 { return gridGeometry_->numDofs(); }
267
269 const Problem& problem() const
270 { return *problem_; }
271
274 { return *gridGeometry_; }
275
277 const GridView& gridView() const
278 { return gridGeometry().gridView(); }
279
282 { return *gridVariables_; }
283
286 { return *gridVariables_; }
287
290 { return *jacobian_; }
291
294 { return *residual_; }
295
297 const SolutionVector& prevSol() const
298 { return *prevSol_; }
299
304 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
305 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
306
312 { prevSol_ = &u; }
313
318 { return isStationaryProblem_; }
319
323 LocalResidual localResidual() const
324 { return LocalResidual(problem_.get(), timeLoop_.get()); }
325
330 {
331 gridVariables().update(cursol);
332 }
333
337 void resetTimeStep(const SolutionVector &cursol)
338 {
339 gridVariables().resetTimeStep(cursol);
340 }
341
342private:
346 void setJacobianPattern_()
347 {
348 // resize the jacobian and the residual
349 const auto numDofs = this->numDofs();
350 jacobian_->setSize(numDofs, numDofs);
351
352 // create occupation pattern of the jacobian
353 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
354
355 // export pattern to jacobian
356 occupationPattern.exportIdx(*jacobian_);
357 }
358
360 void setResidualSize_()
361 { residual_->resize(numDofs()); }
362
364 void maybeComputeColors_()
365 {
366 if (enableMultithreading_)
367 elementSets_ = computeColoring(gridGeometry()).sets;
368 }
369
370 // reset the residual vector to 0.0
371 void resetResidual_()
372 {
373 if(!residual_)
374 {
375 residual_ = std::make_shared<ResidualType>();
376 setResidualSize_();
377 }
378
379 (*residual_) = 0.0;
380 }
381
382 // reset the Jacobian matrix to 0.0
383 template <class PartialReassembler = DefaultPartialReassembler>
384 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
385 {
386 if(!jacobian_)
387 {
388 jacobian_ = std::make_shared<JacobianMatrix>();
389 jacobian_->setBuildMode(JacobianMatrix::random);
390 setJacobianPattern_();
391 }
392
393 if (partialReassembler)
394 partialReassembler->resetJacobian(*this);
395 else
396 *jacobian_ = 0.0;
397 }
398
399 // check if the assembler is in a correct state for assembly
400 void checkAssemblerState_() const
401 {
402 if (!isStationaryProblem_ && !prevSol_)
403 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
404 }
405
411 template<typename AssembleElementFunc>
412 void assemble_(AssembleElementFunc&& assembleElement) const
413 {
414 // a state that will be checked on all processes
415 bool succeeded = false;
416
417 // try assembling using the local assembly function
418 try
419 {
420 if (enableMultithreading_)
421 {
422 assert(elementSets_.size() > 0);
423
424 // make this element loop run in parallel
425 // for this we have to color the elements so that we don't get
426 // race conditions when writing into the global matrix
427 // each color can be assembled using multiple threads
428 for (const auto& elements : elementSets_)
429 {
430 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
431 {
432 const auto element = gridView().grid().entity(elements[i]);
433 assembleElement(element);
434 });
435 }
436 }
437 else
438 for (const auto& element : elements(gridView()))
439 assembleElement(element);
440
441 // if we get here, everything worked well on this process
442 succeeded = true;
443 }
444 // throw exception if a problem occurred
445 catch (NumericalProblem &e)
446 {
447 std::cout << "rank " << gridView().comm().rank()
448 << " caught an exception while assembling:" << e.what()
449 << "\n";
450 succeeded = false;
451 }
452
453 // make sure everything worked well on all processes
454 if (gridView().comm().size() > 1)
455 succeeded = gridView().comm().min(succeeded);
456
457 // if not succeeded rethrow the error on all processes
458 if (!succeeded)
459 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
460 }
461
462 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethods::box, void>
463 enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry)
464 {
465 for (const auto& m : gridGeometry.periodicVertexMap())
466 {
467 if (m.first < m.second)
468 {
469 // add the second row to the first
470 res[m.first] += res[m.second];
471 const auto end = jac[m.second].end();
472 for (auto it = jac[m.second].begin(); it != end; ++it)
473 jac[m.first][it.index()] += (*it);
474
475 // enforce constraint in second row
476 res[m.second] = curSol[m.second] - curSol[m.first];
477 for (auto it = jac[m.second].begin(); it != end; ++it)
478 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
479 }
480 }
481 }
482
483 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethods::box, void>
484 enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry) {}
485
487 std::shared_ptr<const Problem> problem_;
488
490 std::shared_ptr<const GridGeometry> gridGeometry_;
491
493 std::shared_ptr<GridVariables> gridVariables_;
494
496 std::shared_ptr<const TimeLoop> timeLoop_;
497
499 const SolutionVector* prevSol_ = nullptr;
500
502 bool isStationaryProblem_;
503
505 std::shared_ptr<JacobianMatrix> jacobian_;
506 std::shared_ptr<ResidualType> residual_;
507
509 bool enableMultithreading_ = false;
510 std::deque<std::vector<ElementSeed>> elementSets_;
511};
512
513} // namespace Dumux
514
515#endif
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 (cell-centered methods)
Definition assembly/cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition assembly/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:293
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:257
const GridView & gridView() const
The gridview.
Definition assembly/fvassembler.hh:277
const Problem & problem() const
The problem.
Definition assembly/fvassembler.hh:269
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:244
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition assembly/fvassembler.hh:317
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:323
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
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition assembly/fvassembler.hh:265
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:224
GridVariables & gridVariables()
The global grid variables.
Definition assembly/fvassembler.hh:281
JacobianMatrix & jacobian()
The jacobian matrix.
Definition assembly/fvassembler.hh:289
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:304
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:311
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition assembly/fvassembler.hh:273
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition assembly/fvassembler.hh:337
const GridVariables & gridVariables() const
The global grid variables.
Definition assembly/fvassembler.hh:285
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:297
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition assembly/fvassembler.hh:329
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 (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:76
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:57
Traits specifying if a given discretization tag supports coloring.
Definition coloring.hh:292
Type traits to be used with vector types.