3.1-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:
70 using FVGridGeometry [[deprecated("Use GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry;
73
74 using ResidualType = SolutionVector;
75
81 FVAssembler(std::shared_ptr<const Problem> problem,
82 std::shared_ptr<const GridGeometry> gridGeometry,
83 std::shared_ptr<GridVariables> gridVariables)
84 : problem_(problem)
85 , gridGeometry_(gridGeometry)
86 , gridVariables_(gridVariables)
87 , timeLoop_()
88 , isStationaryProblem_(true)
89 {
90 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
91 }
92
97 [[deprecated("Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
98 FVAssembler(std::shared_ptr<const Problem> problem,
99 std::shared_ptr<const GridGeometry> gridGeometry,
100 std::shared_ptr<GridVariables> gridVariables,
101 std::shared_ptr<const TimeLoop> timeLoop)
102 : problem_(problem)
103 , gridGeometry_(gridGeometry)
104 , gridVariables_(gridVariables)
105 , timeLoop_(timeLoop)
106 , isStationaryProblem_(!timeLoop)
107 {}
108
114 FVAssembler(std::shared_ptr<const Problem> problem,
115 std::shared_ptr<const GridGeometry> gridGeometry,
116 std::shared_ptr<GridVariables> gridVariables,
117 std::shared_ptr<const TimeLoop> timeLoop,
118 const SolutionVector& prevSol)
119 : problem_(problem)
120 , gridGeometry_(gridGeometry)
121 , gridVariables_(gridVariables)
122 , timeLoop_(timeLoop)
123 , prevSol_(&prevSol)
124 , isStationaryProblem_(!timeLoop)
125 {}
126
131 template<class PartialReassembler = DefaultPartialReassembler>
132 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
133 {
134 checkAssemblerState_();
135 resetJacobian_(partialReassembler);
136 resetResidual_();
137
138 assemble_([&](const Element& element)
139 {
140 LocalAssembler localAssembler(*this, element, curSol);
141 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
142 });
143
144 enforcePeriodicConstraints_(*jacobian_, *residual_, *gridGeometry_);
145 }
146
150 void assembleJacobian(const SolutionVector& curSol)
151 {
152 checkAssemblerState_();
153 resetJacobian_();
154
155 assemble_([&](const Element& element)
156 {
157 LocalAssembler localAssembler(*this, element, curSol);
158 localAssembler.assembleJacobianAndResidual(*jacobian_, *gridVariables_);
159 });
160 }
161
163 void assembleResidual(const SolutionVector& curSol)
164 {
165 resetResidual_();
166 assembleResidual(*residual_, curSol);
167 }
168
170 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
171 {
172 checkAssemblerState_();
173
174 // update the grid variables for the case of active caching
175 gridVariables_->update(curSol);
176
177 assemble_([&](const Element& element)
178 {
179 LocalAssembler localAssembler(*this, element, curSol);
180 localAssembler.assembleResidual(r);
181 });
182 }
183
185 Scalar residualNorm(const SolutionVector& curSol) const
186 {
188 assembleResidual(residual, curSol);
189
190 // for box communicate the residual with the neighboring processes
191 if (isBox && gridView().comm().size() > 1)
192 {
193 using VertexMapper = typename GridGeometry::VertexMapper;
195 sumResidualHandle(residual, gridGeometry_->vertexMapper());
196 gridView().communicate(sumResidualHandle,
197 Dune::InteriorBorder_InteriorBorder_Interface,
198 Dune::ForwardCommunication);
199 }
200
201 // calculate the square norm of the residual
202 Scalar result2 = residual.two_norm2();
203 if (gridView().comm().size() > 1)
204 result2 = gridView().comm().sum(result2);
205
206 using std::sqrt;
207 return sqrt(result2);
208 }
209
215 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
216 std::shared_ptr<SolutionVector> r)
217 {
218 jacobian_ = A;
219 residual_ = r;
220
221 // check and/or set the BCRS matrix's build mode
222 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
223 jacobian_->setBuildMode(JacobianMatrix::random);
224 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
225 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
226
229 }
230
236 {
237 jacobian_ = std::make_shared<JacobianMatrix>();
238 jacobian_->setBuildMode(JacobianMatrix::random);
239 residual_ = std::make_shared<SolutionVector>();
240
243 }
244
249 {
250 // resize the jacobian and the residual
251 const auto numDofs = this->numDofs();
252 jacobian_->setSize(numDofs, numDofs);
253
254 // create occupation pattern of the jacobian
255 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
256
257 // export pattern to jacobian
258 occupationPattern.exportIdx(*jacobian_);
259 }
260
263 { residual_->resize(numDofs()); }
264
266 std::size_t numDofs() const
267 { return gridGeometry_->numDofs(); }
268
270 const Problem& problem() const
271 { return *problem_; }
272
274 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
276 { return gridGeometry(); }
277
280 { return *gridGeometry_; }
281
283 const GridView& gridView() const
284 { return gridGeometry().gridView(); }
285
288 { return *gridVariables_; }
289
292 { return *gridVariables_; }
293
296 { return *jacobian_; }
297
299 SolutionVector& residual()
300 { return *residual_; }
301
303 const SolutionVector& prevSol() const
304 { return *prevSol_; }
305
310 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
311 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
312
317 void setPreviousSolution(const SolutionVector& u)
318 { prevSol_ = &u; }
319
324 { return isStationaryProblem_; }
325
329 LocalResidual localResidual() const
330 { return LocalResidual(problem_.get(), timeLoop_.get()); }
331
335 void updateGridVariables(const SolutionVector &cursol)
336 {
337 gridVariables().update(cursol);
338 }
339
343 void resetTimeStep(const SolutionVector &cursol)
344 {
345 gridVariables().resetTimeStep(cursol);
346 }
347
348private:
349 // reset the residual vector to 0.0
350 void resetResidual_()
351 {
352 if(!residual_)
353 {
354 residual_ = std::make_shared<SolutionVector>();
356 }
357
358 (*residual_) = 0.0;
359 }
360
361 // reset the Jacobian matrix to 0.0
362 template <class PartialReassembler = DefaultPartialReassembler>
363 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
364 {
365 if(!jacobian_)
366 {
367 jacobian_ = std::make_shared<JacobianMatrix>();
368 jacobian_->setBuildMode(JacobianMatrix::random);
370 }
371
372 if (partialReassembler)
373 partialReassembler->resetJacobian(*this);
374 else
375 *jacobian_ = 0.0;
376 }
377
378 // check if the assembler is in a correct state for assembly
379 void checkAssemblerState_() const
380 {
381 if (!isStationaryProblem_ && !prevSol_)
382 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
383 }
384
390 template<typename AssembleElementFunc>
391 void assemble_(AssembleElementFunc&& assembleElement) const
392 {
393 // a state that will be checked on all processes
394 bool succeeded = false;
395
396 // try assembling using the local assembly function
397 try
398 {
399 // let the local assembler add the element contributions
400 for (const auto& element : elements(gridView()))
401 assembleElement(element);
402
403 // if we get here, everything worked well on this process
404 succeeded = true;
405 }
406 // throw exception if a problem ocurred
407 catch (NumericalProblem &e)
408 {
409 std::cout << "rank " << gridView().comm().rank()
410 << " caught an exception while assembling:" << e.what()
411 << "\n";
412 succeeded = false;
413 }
414
415 // make sure everything worked well on all processes
416 if (gridView().comm().size() > 1)
417 succeeded = gridView().comm().min(succeeded);
418
419 // if not succeeded rethrow the error on all processes
420 if (!succeeded)
421 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
422 }
423
424 template<class GG> std::enable_if_t<GG::discMethod == DiscretizationMethod::box, void>
425 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const GG& gridGeometry)
426 {
427 for (const auto& m : gridGeometry.periodicVertexMap())
428 {
429 if (m.first < m.second)
430 {
431 // add the second row to the first
432 res[m.first] += res[m.second];
433 const auto end = jac[m.second].end();
434 for (auto it = jac[m.second].begin(); it != end; ++it)
435 jac[m.first][it.index()] += (*it);
436
437 // enforce constraint in second row
438 for (auto it = jac[m.second].begin(); it != end; ++it)
439 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
440 }
441 }
442 }
443
444 template<class GG> std::enable_if_t<GG::discMethod != DiscretizationMethod::box, void>
445 enforcePeriodicConstraints_(JacobianMatrix& jac, SolutionVector& res, const GG& gridGeometry) {}
446
448 std::shared_ptr<const Problem> problem_;
449
451 std::shared_ptr<const GridGeometry> gridGeometry_;
452
454 std::shared_ptr<GridVariables> gridVariables_;
455
457 std::shared_ptr<const TimeLoop> timeLoop_;
458
460 const SolutionVector* prevSol_ = nullptr;
461
463 bool isStationaryProblem_;
464
466 std::shared_ptr<JacobianMatrix> jacobian_;
467 std::shared_ptr<SolutionVector> residual_;
468};
469
470} // namespace Dumux
471
472#endif
An assembler for Jacobian and residual contribution per element (box method)
An assembler for Jacobian and residual contribution per element (cell-centered methods)
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.
Manages the handling of time dependent problems.
The available discretization methods in Dumux.
Provides data handles for parallel communication which operate on vertices.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
make the local view function available whenever we use the grid geometry
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
GridGeometry FVGridGeometry
Definition: assembly/fvassembler.hh:70
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:98
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:150
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:283
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:270
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:235
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:74
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:323
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:329
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:114
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:185
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:266
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:287
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:295
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:170
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:71
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:310
void setResidualSize()
Resizes the residual.
Definition: assembly/fvassembler.hh:262
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:163
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:72
void setJacobianPattern()
Resizes the jacobian and sets the jacobian' sparsity pattern.
Definition: assembly/fvassembler.hh:248
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:317
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:279
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:343
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:291
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:215
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:132
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:81
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:299
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:303
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:335
const GridGeometry & fvGridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:275
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
Manages the handling of time dependent problems.
Definition: timeloop.hh:65
Data handle for parallel communication which sums up all values are attached to vertices.
Definition: vertexhandles.hh:42
Declares all properties used in Dumux.