3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/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 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
28
29#include <type_traits>
30#include <tuple>
31
32#include <dune/common/hybridutilities.hh>
33#include <dune/istl/matrixindexset.hh>
34
42
47
49
50namespace Dumux {
51
61template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
63{
64 template<std::size_t id>
65 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
66
67public:
68 using Traits = MDTraits;
69
70 using Scalar = typename MDTraits::Scalar;
71
73 template<std::size_t id>
75
76 template<std::size_t id>
77 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
78
79 template<std::size_t id>
80 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
81
82 template<std::size_t id>
83 using Problem = typename MDTraits::template SubDomain<id>::Problem;
84
85 using JacobianMatrix = typename MDTraits::JacobianMatrix;
86 using SolutionVector = typename MDTraits::SolutionVector;
88
89 using CouplingManager = CMType;
90
94 static constexpr bool isImplicit()
95 { return useImplicitAssembly; }
96
97private:
98
99 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
100 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
101 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
102
104 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
105
106 template<DiscretizationMethod discMethod, std::size_t id>
107 struct SubDomainAssemblerType;
108
109 template<std::size_t id>
110 struct SubDomainAssemblerType<DiscretizationMethod::cctpfa, id>
111 {
112 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
113 };
114
115 template<std::size_t id>
116 struct SubDomainAssemblerType<DiscretizationMethod::ccmpfa, id>
117 {
118 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
119 };
120
121 template<std::size_t id>
122 struct SubDomainAssemblerType<DiscretizationMethod::box, id>
123 {
124 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
125 };
126
127 template<std::size_t id>
128 struct SubDomainAssemblerType<DiscretizationMethod::staggered, id>
129 {
130 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
131 };
132
133 template<std::size_t id>
134 using SubDomainAssembler = typename SubDomainAssemblerType<GridGeometry<id>::discMethod, id>::type;
135
136public:
137
138
145 GridGeometryTuple&& gridGeometry,
146 GridVariablesTuple&& gridVariables,
147 std::shared_ptr<CouplingManager> couplingManager)
149 , problemTuple_(problem)
150 , gridGeometryTuple_(gridGeometry)
151 , gridVariablesTuple_(gridVariables)
152 , timeLoop_()
153 , isStationaryProblem_(true)
154 , warningIssued_(false)
155 {
156 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
157 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
158 }
159
166 GridGeometryTuple&& gridGeometry,
167 GridVariablesTuple&& gridVariables,
168 std::shared_ptr<CouplingManager> couplingManager,
169 std::shared_ptr<const TimeLoop> timeLoop,
170 const SolutionVector& prevSol)
172 , problemTuple_(problem)
173 , gridGeometryTuple_(gridGeometry)
174 , gridVariablesTuple_(gridVariables)
175 , timeLoop_(timeLoop)
176 , prevSol_(&prevSol)
177 , isStationaryProblem_(false)
178 , warningIssued_(false)
179 {
180 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
181 }
182
188 {
189 checkAssemblerState_();
190 resetJacobian_();
191 resetResidual_();
192
193 using namespace Dune::Hybrid;
194 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
195 {
196 auto& jacRow = (*jacobian_)[domainId];
197 auto& subRes = (*residual_)[domainId];
198 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
199 });
200 }
201
204 {
205 resetResidual_();
206 assembleResidual(*residual_, curSol);
207 }
208
211 {
212 r = 0.0;
213
214 checkAssemblerState_();
215
216 // update the grid variables for the case of active caching
217 updateGridVariables(curSol);
218
219 using namespace Dune::Hybrid;
220 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
221 {
222 auto& subRes = r[domainId];
223 this->assembleResidual_(domainId, subRes, curSol);
224 });
225 }
226
229 {
232 assembleResidual(residual, curSol);
233
234 // calculate the squared norm of the residual
235 Scalar resultSquared = 0.0;
236
237 // for box communicate the residual with the neighboring processes
238 using namespace Dune::Hybrid;
239 forEach(integralRange(Dune::Hybrid::size(residual)), [&](const auto domainId)
240 {
241 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
242 const auto& gridView = gridGeometry->gridView();
243
244 if (gridView.comm().size() > 1 && gridView.overlapSize(0) == 0)
245 {
247 {
248 using GV = typename GridGeometry<domainId>::GridView;
249 using DM = typename GridGeometry<domainId>::VertexMapper;
251
252 PVHelper vectorHelper(gridView, gridGeometry->vertexMapper());
253
254 vectorHelper.makeNonOverlappingConsistent(residual[domainId]);
255 }
256 }
257 else if (!warningIssued_)
258 {
259 if (gridView.comm().size() > 1 && gridView.comm().rank() == 0)
260 std::cout << "\nWarning: norm calculation adds entries corresponding to\n"
261 << "overlapping entities multiple times. Please use the norm\n"
262 << "function provided by a linear solver instead." << std::endl;
263
264 warningIssued_ = true;
265 }
266
267 Scalar localNormSquared = residual[domainId].two_norm2();
268
269 if (gridView.comm().size() > 1)
270 {
271 localNormSquared = gridView.comm().sum(localNormSquared);
272 }
273
274 resultSquared += localNormSquared;
275 });
276
277 using std::sqrt;
278 return sqrt(resultSquared);
279 }
280
286 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
287 std::shared_ptr<SolutionVector> r)
288 {
289 jacobian_ = A;
290 residual_ = r;
291
292 setJacobianBuildMode(*jacobian_);
293 setJacobianPattern(*jacobian_);
294 setResidualSize(*residual_);
295 }
296
302 {
303 jacobian_ = std::make_shared<JacobianMatrix>();
304 residual_ = std::make_shared<SolutionVector>();
305
306 setJacobianBuildMode(*jacobian_);
307 setJacobianPattern(*jacobian_);
308 setResidualSize(*residual_);
309 }
310
315 {
316 using namespace Dune::Hybrid;
317 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
318 {
319 forEach(jac[i], [&](auto& jacBlock)
320 {
321 using BlockType = std::decay_t<decltype(jacBlock)>;
322 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
323 jacBlock.setBuildMode(BlockType::BuildMode::random);
324 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
325 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
326 });
327 });
328 }
329
334 {
335 using namespace Dune::Hybrid;
336 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
337 {
338 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
339 {
340 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
341 pattern.exportIdx(jac[domainI][domainJ]);
342 });
343 });
344 }
345
350 {
351 using namespace Dune::Hybrid;
352 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
353 { res[domainId].resize(this->numDofs(domainId)); });
354 }
355
360 {
361 using namespace Dune::Hybrid;
362 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
363 { this->gridVariables(domainId).update(curSol[domainId]); });
364 }
365
369 void resetTimeStep(const SolutionVector& curSol)
370 {
371 using namespace Dune::Hybrid;
372 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
373 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
374 }
375
377 template<std::size_t i>
378 std::size_t numDofs(Dune::index_constant<i> domainId) const
379 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
380
382 template<std::size_t i>
383 const auto& problem(Dune::index_constant<i> domainId) const
384 { return *std::get<domainId>(problemTuple_); }
385
387 template<std::size_t i>
388 const auto& gridGeometry(Dune::index_constant<i> domainId) const
389 { return *std::get<domainId>(gridGeometryTuple_); }
390
392 template<std::size_t i>
393 const auto& gridView(Dune::index_constant<i> domainId) const
394 { return gridGeometry(domainId).gridView(); }
395
397 template<std::size_t i>
398 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
399 { return *std::get<domainId>(gridVariablesTuple_); }
400
402 template<std::size_t i>
403 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
404 { return *std::get<domainId>(gridVariablesTuple_); }
405
408 { return *couplingManager_; }
409
412 { return *jacobian_; }
413
416 { return *residual_; }
417
419 const SolutionVector& prevSol() const
420 { return *prevSol_; }
421
426 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
427 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
428
434 { prevSol_ = &u; }
435
440 { return isStationaryProblem_; }
441
445 template<std::size_t i>
446 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
447 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
448
449protected:
451 std::shared_ptr<CouplingManager> couplingManager_;
452
453private:
454 // reset the residual vector to 0.0
455 void resetResidual_()
456 {
457 if(!residual_)
458 {
459 residual_ = std::make_shared<SolutionVector>();
460 setResidualSize(*residual_);
461 }
462
463 (*residual_) = 0.0;
464 }
465
466 // reset the jacobian vector to 0.0
467 void resetJacobian_()
468 {
469 if(!jacobian_)
470 {
471 jacobian_ = std::make_shared<JacobianMatrix>();
472 setJacobianBuildMode(*jacobian_);
473 setJacobianPattern(*jacobian_);
474 }
475
476 (*jacobian_) = 0.0;
477 }
478
479 // check if the assembler is in a correct state for assembly
480 void checkAssemblerState_() const
481 {
482 if (!isStationaryProblem_ && !prevSol_)
483 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
484
485 if (isStationaryProblem_ && prevSol_)
486 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
487 << " Did you forget to set the timeLoop to make this problem instationary?");
488 }
489
490 template<std::size_t i, class JacRow, class SubRes>
491 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
492 const SolutionVector& curSol)
493 {
494 assemble_(domainId, [&](const auto& element)
495 {
496 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
497 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
498 });
499 }
500
501 template<std::size_t i, class SubRes>
502 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
503 const SolutionVector& curSol)
504 {
505 assemble_(domainId, [&](const auto& element)
506 {
507 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
508 subDomainAssembler.assembleResidual(subRes);
509 });
510 }
511
518 template<std::size_t i, class AssembleElementFunc>
519 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
520 {
521 // let the local assembler add the element contributions
522 for (const auto& element : elements(gridView(domainId)))
523 assembleElement(element);
524 }
525
526 // get diagonal block pattern
527 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
528 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
529 Dune::index_constant<j> domainJ) const
530 {
531 const auto& gg = gridGeometry(domainI);
532 auto pattern = getJacobianPattern<isImplicit()>(gg);
533 couplingManager_->extendJacobianPattern(domainI, pattern);
534 return pattern;
535 }
536
537 // get coupling block pattern
538 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
539 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
540 Dune::index_constant<j> domainJ) const
541 {
542 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
543 domainI, gridGeometry(domainI),
544 domainJ, gridGeometry(domainJ));
545 }
546
548 ProblemTuple problemTuple_;
549
551 GridGeometryTuple gridGeometryTuple_;
552
554 GridVariablesTuple gridVariablesTuple_;
555
557 std::shared_ptr<const TimeLoop> timeLoop_;
558
560 const SolutionVector* prevSol_ = nullptr;
561
563 bool isStationaryProblem_;
564
566 std::shared_ptr<JacobianMatrix> jacobian_;
567 std::shared_ptr<SolutionVector> residual_;
568
570 bool warningIssued_;
571};
572
573} // end namespace Dumux
574
575#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.
Utilities for template meta programming.
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
Helper function to generate Jacobian pattern for multi domain models.
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
An assembler for Jacobian and residual contribution per element (box methods) for multidomain problem...
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
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
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Definition: common/properties.hh:74
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
Definition: parallelhelpers.hh:430
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:63
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:378
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:359
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:86
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:87
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:83
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:411
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:94
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:210
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:85
CMType CouplingManager
Definition: multidomain/fvassembler.hh:89
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:203
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:369
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:228
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:403
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: multidomain/fvassembler.hh:433
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:77
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: multidomain/fvassembler.hh:286
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:388
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:349
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:301
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:393
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: multidomain/fvassembler.hh:165
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:80
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:407
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:426
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:419
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:398
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:451
MDTraits Traits
Definition: multidomain/fvassembler.hh:68
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:314
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:74
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:144
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:446
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:383
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:187
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:439
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:70
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:333
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:415
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Declares all properties used in Dumux.
Manages the handling of time dependent problems.