Loading [MathJax]/extensions/tex2jax.js
3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 // issue a warning if the caluclation is used in parallel with overlap
171 static bool warningIssued = false;
172
173 if (gridView().comm().size() > 1 && gridView().overlapSize(0) == 0)
174 {
175 if constexpr (isBox)
176 {
177 using DM = typename GridGeometry::VertexMapper;
179
180 PVHelper vectorHelper(gridView(), gridGeometry_->vertexMapper());
181
182 vectorHelper.makeNonOverlappingConsistent(residual);
183 }
184 }
185 else if (!warningIssued)
186 {
187 if (gridView().comm().size() > 1 && gridView().comm().rank() == 0)
188 std::cout << "\nWarning: norm calculation adds entries corresponding to\n"
189 << "overlapping entities multiple times. Please use the norm\n"
190 << "function provided by a linear solver instead." << std::endl;
191
192 warningIssued = true;
193 }
194
195 // calculate the square norm of the residual
196 Scalar result2 = residual.two_norm2();
197 if (gridView().comm().size() > 1)
198 result2 = gridView().comm().sum(result2);
199
200 using std::sqrt;
201 return sqrt(result2);
202 }
203
209 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
210 std::shared_ptr<SolutionVector> r)
211 {
212 jacobian_ = A;
213 residual_ = r;
214
215 // check and/or set the BCRS matrix's build mode
216 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
217 jacobian_->setBuildMode(JacobianMatrix::random);
218 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
219 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
220
223 }
224
230 {
231 jacobian_ = std::make_shared<JacobianMatrix>();
232 jacobian_->setBuildMode(JacobianMatrix::random);
233 residual_ = std::make_shared<SolutionVector>();
234
237 }
238
243 {
244 // resize the jacobian and the residual
245 const auto numDofs = this->numDofs();
246 jacobian_->setSize(numDofs, numDofs);
247
248 // create occupation pattern of the jacobian
249 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
250
251 // export pattern to jacobian
252 occupationPattern.exportIdx(*jacobian_);
253 }
254
257 { residual_->resize(numDofs()); }
258
260 std::size_t numDofs() const
261 { return gridGeometry_->numDofs(); }
262
264 const Problem& problem() const
265 { return *problem_; }
266
269 { return *gridGeometry_; }
270
272 const GridView& gridView() const
273 { return gridGeometry().gridView(); }
274
277 { return *gridVariables_; }
278
281 { return *gridVariables_; }
282
285 { return *jacobian_; }
286
288 SolutionVector& residual()
289 { return *residual_; }
290
292 const SolutionVector& prevSol() const
293 { return *prevSol_; }
294
299 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
300 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
301
306 void setPreviousSolution(const SolutionVector& u)
307 { prevSol_ = &u; }
308
313 { return isStationaryProblem_; }
314
318 LocalResidual localResidual() const
319 { return LocalResidual(problem_.get(), timeLoop_.get()); }
320
324 void updateGridVariables(const SolutionVector &cursol)
325 {
326 gridVariables().update(cursol);
327 }
328
332 void resetTimeStep(const SolutionVector &cursol)
333 {
334 gridVariables().resetTimeStep(cursol);
335 }
336
337private:
338 // reset the residual vector to 0.0
339 void resetResidual_()
340 {
341 if(!residual_)
342 {
343 residual_ = std::make_shared<SolutionVector>();
345 }
346
347 (*residual_) = 0.0;
348 }
349
350 // reset the Jacobian matrix to 0.0
351 template <class PartialReassembler = DefaultPartialReassembler>
352 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
353 {
354 if(!jacobian_)
355 {
356 jacobian_ = std::make_shared<JacobianMatrix>();
357 jacobian_->setBuildMode(JacobianMatrix::random);
359 }
360
361 if (partialReassembler)
362 partialReassembler->resetJacobian(*this);
363 else
364 *jacobian_ = 0.0;
365 }
366
367 // check if the assembler is in a correct state for assembly
368 void checkAssemblerState_() const
369 {
370 if (!isStationaryProblem_ && !prevSol_)
371 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
372 }
373
379 template<typename AssembleElementFunc>
380 void assemble_(AssembleElementFunc&& assembleElement) const
381 {
382 // a state that will be checked on all processes
383 bool succeeded = false;
384
385 // try assembling using the local assembly function
386 try
387 {
388 // let the local assembler add the element contributions
389 for (const auto& element : elements(gridView()))
390 assembleElement(element);
391
392 // if we get here, everything worked well on this process
393 succeeded = true;
394 }
395 // throw exception if a problem ocurred
396 catch (NumericalProblem &e)
397 {
398 std::cout << "rank " << gridView().comm().rank()
399 << " caught an exception while assembling:" << e.what()
400 << "\n";
401 succeeded = false;
402 }
403
404 // make sure everything worked well on all processes
405 if (gridView().comm().size() > 1)
406 succeeded = gridView().comm().min(succeeded);
407
408 // if not succeeded rethrow the error on all processes
409 if (!succeeded)
410 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
411 }
412
413 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
414 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry)
415 {
416 for (const auto& m : gridGeometry.periodicVertexMap())
417 {
418 if (m.first < m.second)
419 {
420 // add the second row to the first
421 res[m.first] += res[m.second];
422 const auto end = jac[m.second].end();
423 for (auto it = jac[m.second].begin(); it != end; ++it)
424 jac[m.first][it.index()] += (*it);
425
426 // enforce constraint in second row
427 res[m.second] = curSol[m.second] - curSol[m.first];
428 for (auto it = jac[m.second].begin(); it != end; ++it)
429 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
430 }
431 }
432 }
433
434 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
435 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const SolutionVector& curSol, const GG& gridGeometry) {}
436
438 std::shared_ptr<const Problem> problem_;
439
441 std::shared_ptr<const GridGeometry> gridGeometry_;
442
444 std::shared_ptr<GridVariables> gridVariables_;
445
447 std::shared_ptr<const TimeLoop> timeLoop_;
448
450 const SolutionVector* prevSol_ = nullptr;
451
453 bool isStationaryProblem_;
454
456 std::shared_ptr<JacobianMatrix> jacobian_;
457 std::shared_ptr<SolutionVector> residual_;
458};
459
460} // namespace Dumux
461
462#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)
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
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
Definition: propertysystem.hh:150
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: cclocalassembler.hh:136
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:272
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:264
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:229
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:73
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:312
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:318
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:260
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:276
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:284
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:299
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:256
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:242
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:306
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:268
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:332
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:280
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:209
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:288
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:292
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:324
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
Definition: parallelhelpers.hh:430
Declares all properties used in Dumux.
Manages the handling of time dependent problems.