3.3.0
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 <type_traits>
28
29#include <dune/istl/matrixindexset.hh>
30
35
36#include "jacobianpattern.hh"
37#include "diffmethod.hh"
38#include "boxlocalassembler.hh"
39#include "cclocalassembler.hh"
40
41namespace Dumux {
42
50template<class TypeTag, DiffMethod diffMethod, bool isImplicit = true>
52{
55 using Element = typename GridView::template Codim<0>::Entity;
58
60 static constexpr bool isBox = discMethod == DiscretizationMethod::box;
61
63 using LocalAssembler = std::conditional_t<isBox, BoxLocalAssembler<TypeTag, ThisType, diffMethod, isImplicit>,
65
66public:
72
73 using ResidualType = SolutionVector;
74
80 FVAssembler(std::shared_ptr<const Problem> problem,
81 std::shared_ptr<const GridGeometry> gridGeometry,
82 std::shared_ptr<GridVariables> gridVariables)
83 : problem_(problem)
84 , gridGeometry_(gridGeometry)
85 , gridVariables_(gridVariables)
86 , timeLoop_()
87 , isStationaryProblem_(true)
88 {
89 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
90 }
91
97 FVAssembler(std::shared_ptr<const Problem> problem,
98 std::shared_ptr<const GridGeometry> gridGeometry,
99 std::shared_ptr<GridVariables> gridVariables,
100 std::shared_ptr<const TimeLoop> timeLoop,
101 const SolutionVector& prevSol)
102 : problem_(problem)
103 , gridGeometry_(gridGeometry)
104 , gridVariables_(gridVariables)
105 , timeLoop_(timeLoop)
106 , prevSol_(&prevSol)
107 , isStationaryProblem_(!timeLoop)
108 {}
109
114 template<class PartialReassembler = DefaultPartialReassembler>
115 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
116 {
117 checkAssemblerState_();
118 resetJacobian_(partialReassembler);
119 resetResidual_();
120
121 assemble_([&](const Element& element)
122 {
123 LocalAssembler localAssembler(*this, element, curSol);
124 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
125 });
126
127 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
128 }
129
133 void assembleJacobian(const SolutionVector& curSol)
134 {
135 checkAssemblerState_();
136 resetJacobian_();
137
138 assemble_([&](const Element& element)
139 {
140 LocalAssembler localAssembler(*this, element, curSol);
141 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
142 });
143 }
144
146 void assembleResidual(const SolutionVector& curSol)
147 {
148 resetResidual_();
149 assembleResidual(*residual_, curSol);
150 }
151
153 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
154 {
155 checkAssemblerState_();
156
157 assemble_([&](const Element& element)
158 {
159 LocalAssembler localAssembler(*this, element, curSol);
160 localAssembler.assembleResidual(r);
161 });
162 }
163
165 Scalar residualNorm(const SolutionVector& curSol) const
166 {
168 assembleResidual(residual, curSol);
169
170 // for box communicate the residual with the neighboring processes
171 if (isBox && gridView().comm().size() > 1)
172 {
173 using VertexMapper = typename GridGeometry::VertexMapper;
175 sumResidualHandle(gridGeometry_->vertexMapper(), residual);
176 gridView().communicate(sumResidualHandle,
177 Dune::InteriorBorder_InteriorBorder_Interface,
178 Dune::ForwardCommunication);
179 }
180
181 // calculate the square norm of the residual
182 Scalar result2 = residual.two_norm2();
183 if (gridView().comm().size() > 1)
184 result2 = gridView().comm().sum(result2);
185
186 using std::sqrt;
187 return sqrt(result2);
188 }
189
195 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
196 std::shared_ptr<SolutionVector> r)
197 {
198 jacobian_ = A;
199 residual_ = r;
200
201 // check and/or set the BCRS matrix's build mode
202 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
203 jacobian_->setBuildMode(JacobianMatrix::random);
204 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
205 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
206
209 }
210
216 {
217 jacobian_ = std::make_shared<JacobianMatrix>();
218 jacobian_->setBuildMode(JacobianMatrix::random);
219 residual_ = std::make_shared<SolutionVector>();
220
223 }
224
229 {
230 // resize the jacobian and the residual
231 const auto numDofs = this->numDofs();
232 jacobian_->setSize(numDofs, numDofs);
233
234 // create occupation pattern of the jacobian
235 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
236
237 // export pattern to jacobian
238 occupationPattern.exportIdx(*jacobian_);
239 }
240
243 { residual_->resize(numDofs()); }
244
246 std::size_t numDofs() const
247 { return gridGeometry_->numDofs(); }
248
250 const Problem& problem() const
251 { return *problem_; }
252
255 { return *gridGeometry_; }
256
258 const GridView& gridView() const
259 { return gridGeometry().gridView(); }
260
263 { return *gridVariables_; }
264
267 { return *gridVariables_; }
268
271 { return *jacobian_; }
272
274 SolutionVector& residual()
275 { return *residual_; }
276
278 const SolutionVector& prevSol() const
279 { return *prevSol_; }
280
285 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
286 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
287
292 void setPreviousSolution(const SolutionVector& u)
293 { prevSol_ = &u; }
294
299 { return isStationaryProblem_; }
300
304 LocalResidual localResidual() const
305 { return LocalResidual(problem_.get(), timeLoop_.get()); }
306
310 void updateGridVariables(const SolutionVector &cursol)
311 {
312 gridVariables().update(cursol);
313 }
314
318 void resetTimeStep(const SolutionVector &cursol)
319 {
320 gridVariables().resetTimeStep(cursol);
321 }
322
323private:
324 // reset the residual vector to 0.0
325 void resetResidual_()
326 {
327 if(!residual_)
328 {
329 residual_ = std::make_shared<SolutionVector>();
331 }
332
333 (*residual_) = 0.0;
334 }
335
336 // reset the Jacobian matrix to 0.0
337 template <class PartialReassembler = DefaultPartialReassembler>
338 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
339 {
340 if(!jacobian_)
341 {
342 jacobian_ = std::make_shared<JacobianMatrix>();
343 jacobian_->setBuildMode(JacobianMatrix::random);
345 }
346
347 if (partialReassembler)
348 partialReassembler->resetJacobian(*this);
349 else
350 *jacobian_ = 0.0;
351 }
352
353 // check if the assembler is in a correct state for assembly
354 void checkAssemblerState_() const
355 {
356 if (!isStationaryProblem_ && !prevSol_)
357 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
358 }
359
365 template<typename AssembleElementFunc>
366 void assemble_(AssembleElementFunc&& assembleElement) const
367 {
368 // a state that will be checked on all processes
369 bool succeeded = false;
370
371 // try assembling using the local assembly function
372 try
373 {
374 // let the local assembler add the element contributions
375 for (const auto& element : elements(gridView()))
376 assembleElement(element);
377
378 // if we get here, everything worked well on this process
379 succeeded = true;
380 }
381 // throw exception if a problem ocurred
382 catch (NumericalProblem &e)
383 {
384 std::cout << "rank " << gridView().comm().rank()
385 << " caught an exception while assembling:" << e.what()
386 << "\n";
387 succeeded = false;
388 }
389
390 // make sure everything worked well on all processes
391 if (gridView().comm().size() > 1)
392 succeeded = gridView().comm().min(succeeded);
393
394 // if not succeeded rethrow the error on all processes
395 if (!succeeded)
396 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
397 }
398
399 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
400 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry)
401 {
402 for (const auto& m : gridGeometry.periodicVertexMap())
403 {
404 if (m.first < m.second)
405 {
406 // add the second row to the first
407 res[m.first] += res[m.second];
408 const auto end = jac[m.second].end();
409 for (auto it = jac[m.second].begin(); it != end; ++it)
410 jac[m.first][it.index()] += (*it);
411
412 // enforce constraint in second row
413 res[m.second] = curSol[m.second] - curSol[m.first];
414 for (auto it = jac[m.second].begin(); it != end; ++it)
415 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
416 }
417 }
418 }
419
420 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
421 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry) {}
422
424 std::shared_ptr<const Problem> problem_;
425
427 std::shared_ptr<const GridGeometry> gridGeometry_;
428
430 std::shared_ptr<GridVariables> gridVariables_;
431
433 std::shared_ptr<const TimeLoop> timeLoop_;
434
436 const SolutionVector* prevSol_ = nullptr;
437
439 bool isStationaryProblem_;
440
442 std::shared_ptr<JacobianMatrix> jacobian_;
443 std::shared_ptr<SolutionVector> residual_;
444};
445
446} // namespace Dumux
447
448#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper function to generate Jacobian pattern for different discretization methods.
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
The available discretization methods in Dumux.
Contains a class to exchange entries of a vector.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:121
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:52
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:67
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:133
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:258
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:250
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:215
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:73
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:298
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:68
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:304
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:97
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: assembly/fvassembler.hh:69
Scalar residualNorm(const SolutionVector &curSol) const
compute the residual and return it's vector norm
Definition: assembly/fvassembler.hh:165
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:246
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:262
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:270
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:153
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:70
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:285
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:242
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:146
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:71
void setJacobianPattern()
Resizes the jacobian and sets the jacobian' sparsity pattern.
Definition: assembly/fvassembler.hh:228
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:292
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:254
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:318
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:266
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:195
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: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:80
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:274
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:278
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:310
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:425
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78
Declares all properties used in Dumux.
Manages the handling of time dependent problems.