3.2-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 <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
96 [[deprecated("Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
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 : problem_(problem)
102 , gridGeometry_(gridGeometry)
103 , gridVariables_(gridVariables)
104 , timeLoop_(timeLoop)
105 , isStationaryProblem_(!timeLoop)
106 {}
107
113 FVAssembler(std::shared_ptr<const Problem> problem,
114 std::shared_ptr<const GridGeometry> gridGeometry,
115 std::shared_ptr<GridVariables> gridVariables,
116 std::shared_ptr<const TimeLoop> timeLoop,
117 const SolutionVector& prevSol)
118 : problem_(problem)
119 , gridGeometry_(gridGeometry)
120 , gridVariables_(gridVariables)
121 , timeLoop_(timeLoop)
122 , prevSol_(&prevSol)
123 , isStationaryProblem_(!timeLoop)
124 {}
125
130 template<class PartialReassembler = DefaultPartialReassembler>
131 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
132 {
133 checkAssemblerState_();
134 resetJacobian_(partialReassembler);
135 resetResidual_();
136
137 assemble_([&](const Element& element)
138 {
139 LocalAssembler localAssembler(*this, element, curSol);
140 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
141 });
142
143 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
144 }
145
149 void assembleJacobian(const SolutionVector& curSol)
150 {
151 checkAssemblerState_();
152 resetJacobian_();
153
154 assemble_([&](const Element& element)
155 {
156 LocalAssembler localAssembler(*this, element, curSol);
157 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
158 });
159 }
160
162 void assembleResidual(const SolutionVector& curSol)
163 {
164 resetResidual_();
165 assembleResidual(*residual_, curSol);
166 }
167
169 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
170 {
171 checkAssemblerState_();
172
173 // update the grid variables for the case of active caching
174 gridVariables_->update(curSol);
175
176 assemble_([&](const Element& element)
177 {
178 LocalAssembler localAssembler(*this, element, curSol);
179 localAssembler.assembleResidual(r);
180 });
181 }
182
184 Scalar residualNorm(const SolutionVector& curSol) const
185 {
187 assembleResidual(residual, curSol);
188
189 // for box communicate the residual with the neighboring processes
190 if (isBox && gridView().comm().size() > 1)
191 {
192 using VertexMapper = typename GridGeometry::VertexMapper;
194 sumResidualHandle(gridGeometry_->vertexMapper(), residual);
195 gridView().communicate(sumResidualHandle,
196 Dune::InteriorBorder_InteriorBorder_Interface,
197 Dune::ForwardCommunication);
198 }
199
200 // calculate the square norm of the residual
201 Scalar result2 = residual.two_norm2();
202 if (gridView().comm().size() > 1)
203 result2 = gridView().comm().sum(result2);
204
205 using std::sqrt;
206 return sqrt(result2);
207 }
208
214 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
215 std::shared_ptr<SolutionVector> r)
216 {
217 jacobian_ = A;
218 residual_ = r;
219
220 // check and/or set the BCRS matrix's build mode
221 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
222 jacobian_->setBuildMode(JacobianMatrix::random);
223 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
224 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
225
228 }
229
235 {
236 jacobian_ = std::make_shared<JacobianMatrix>();
237 jacobian_->setBuildMode(JacobianMatrix::random);
238 residual_ = std::make_shared<SolutionVector>();
239
242 }
243
248 {
249 // resize the jacobian and the residual
250 const auto numDofs = this->numDofs();
251 jacobian_->setSize(numDofs, numDofs);
252
253 // create occupation pattern of the jacobian
254 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
255
256 // export pattern to jacobian
257 occupationPattern.exportIdx(*jacobian_);
258 }
259
262 { residual_->resize(numDofs()); }
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
293 SolutionVector& residual()
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
311 void setPreviousSolution(const SolutionVector& u)
312 { prevSol_ = &u; }
313
318 { return isStationaryProblem_; }
319
323 LocalResidual localResidual() const
324 { return LocalResidual(problem_.get(), timeLoop_.get()); }
325
329 void updateGridVariables(const SolutionVector &cursol)
330 {
331 gridVariables().update(cursol);
332 }
333
337 void resetTimeStep(const SolutionVector &cursol)
338 {
339 gridVariables().resetTimeStep(cursol);
340 }
341
342private:
343 // reset the residual vector to 0.0
344 void resetResidual_()
345 {
346 if(!residual_)
347 {
348 residual_ = std::make_shared<SolutionVector>();
350 }
351
352 (*residual_) = 0.0;
353 }
354
355 // reset the Jacobian matrix to 0.0
356 template <class PartialReassembler = DefaultPartialReassembler>
357 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
358 {
359 if(!jacobian_)
360 {
361 jacobian_ = std::make_shared<JacobianMatrix>();
362 jacobian_->setBuildMode(JacobianMatrix::random);
364 }
365
366 if (partialReassembler)
367 partialReassembler->resetJacobian(*this);
368 else
369 *jacobian_ = 0.0;
370 }
371
372 // check if the assembler is in a correct state for assembly
373 void checkAssemblerState_() const
374 {
375 if (!isStationaryProblem_ && !prevSol_)
376 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
377 }
378
384 template<typename AssembleElementFunc>
385 void assemble_(AssembleElementFunc&& assembleElement) const
386 {
387 // a state that will be checked on all processes
388 bool succeeded = false;
389
390 // try assembling using the local assembly function
391 try
392 {
393 // let the local assembler add the element contributions
394 for (const auto& element : elements(gridView()))
395 assembleElement(element);
396
397 // if we get here, everything worked well on this process
398 succeeded = true;
399 }
400 // throw exception if a problem ocurred
401 catch (NumericalProblem &e)
402 {
403 std::cout << "rank " << gridView().comm().rank()
404 << " caught an exception while assembling:" << e.what()
405 << "\n";
406 succeeded = false;
407 }
408
409 // make sure everything worked well on all processes
410 if (gridView().comm().size() > 1)
411 succeeded = gridView().comm().min(succeeded);
412
413 // if not succeeded rethrow the error on all processes
414 if (!succeeded)
415 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
416 }
417
418 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
419 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry)
420 {
421 for (const auto& m : gridGeometry.periodicVertexMap())
422 {
423 if (m.first < m.second)
424 {
425 // add the second row to the first
426 res[m.first] += res[m.second];
427 const auto end = jac[m.second].end();
428 for (auto it = jac[m.second].begin(); it != end; ++it)
429 jac[m.first][it.index()] += (*it);
430
431 // enforce constraint in second row
432 res[m.second] = curSol[m.second] - curSol[m.first];
433 for (auto it = jac[m.second].begin(); it != end; ++it)
434 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
435 }
436 }
437 }
438
439 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
440 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry) {}
441
443 std::shared_ptr<const Problem> problem_;
444
446 std::shared_ptr<const GridGeometry> gridGeometry_;
447
449 std::shared_ptr<GridVariables> gridVariables_;
450
452 std::shared_ptr<const TimeLoop> timeLoop_;
453
455 const SolutionVector* prevSol_ = nullptr;
456
458 bool isStationaryProblem_;
459
461 std::shared_ptr<JacobianMatrix> jacobian_;
462 std::shared_ptr<SolutionVector> residual_;
463};
464
465} // namespace Dumux
466
467#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)
Manages the handling of time dependent problems.
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:123
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:52
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop)
The constructor for instationary problems.
Definition: assembly/fvassembler.hh:97
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:149
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:277
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:269
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:234
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:73
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:68
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:323
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:113
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:184
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:265
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:169
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:304
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:261
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:162
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:247
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 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:214
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:131
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:293
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
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:79
Declares all properties used in Dumux.