3.3.0
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
41
46
48
49namespace Dumux {
50
60template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
62{
63 template<std::size_t id>
64 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
65
66public:
67 using Traits = MDTraits;
68
69 using Scalar = typename MDTraits::Scalar;
70
72 template<std::size_t id>
74
75 template<std::size_t id>
76 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
77
78 template<std::size_t id>
79 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
80
81 template<std::size_t id>
82 using Problem = typename MDTraits::template SubDomain<id>::Problem;
83
84 using JacobianMatrix = typename MDTraits::JacobianMatrix;
85 using SolutionVector = typename MDTraits::SolutionVector;
87
88 using CouplingManager = CMType;
89
93 static constexpr bool isImplicit()
94 { return useImplicitAssembly; }
95
96private:
97
98 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
99 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
100 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
101
103 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
104
105 template<DiscretizationMethod discMethod, std::size_t id>
106 struct SubDomainAssemblerType;
107
108 template<std::size_t id>
109 struct SubDomainAssemblerType<DiscretizationMethod::cctpfa, id>
110 {
111 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
112 };
113
114 template<std::size_t id>
115 struct SubDomainAssemblerType<DiscretizationMethod::ccmpfa, id>
116 {
117 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
118 };
119
120 template<std::size_t id>
121 struct SubDomainAssemblerType<DiscretizationMethod::box, id>
122 {
123 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
124 };
125
126 template<std::size_t id>
127 struct SubDomainAssemblerType<DiscretizationMethod::staggered, id>
128 {
129 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
130 };
131
132 template<std::size_t id>
133 using SubDomainAssembler = typename SubDomainAssemblerType<GridGeometry<id>::discMethod, id>::type;
134
135public:
136
137
144 GridGeometryTuple&& gridGeometry,
145 GridVariablesTuple&& gridVariables,
146 std::shared_ptr<CouplingManager> couplingManager)
148 , problemTuple_(problem)
149 , gridGeometryTuple_(gridGeometry)
150 , gridVariablesTuple_(gridVariables)
151 , timeLoop_()
152 , isStationaryProblem_(true)
153 {
154 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
155 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
156 }
157
164 GridGeometryTuple&& gridGeometry,
165 GridVariablesTuple&& gridVariables,
166 std::shared_ptr<CouplingManager> couplingManager,
167 std::shared_ptr<const TimeLoop> timeLoop,
168 const SolutionVector& prevSol)
170 , problemTuple_(problem)
171 , gridGeometryTuple_(gridGeometry)
172 , gridVariablesTuple_(gridVariables)
173 , timeLoop_(timeLoop)
174 , prevSol_(&prevSol)
175 , isStationaryProblem_(false)
176 {
177 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
178 }
179
185 {
186 checkAssemblerState_();
187 resetJacobian_();
188 resetResidual_();
189
190 using namespace Dune::Hybrid;
191 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
192 {
193 auto& jacRow = (*jacobian_)[domainId];
194 auto& subRes = (*residual_)[domainId];
195 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
196 });
197 }
198
201 {
202 resetResidual_();
203 assembleResidual(*residual_, curSol);
204 }
205
208 {
209 checkAssemblerState_();
210
211 // update the grid variables for the case of active caching
212 updateGridVariables(curSol);
213
214 using namespace Dune::Hybrid;
215 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
216 {
217 auto& subRes = r[domainId];
218 this->assembleResidual_(domainId, subRes, curSol);
219 });
220 }
221
225 {
228 assembleResidual(residual, curSol);
229
230 // calculate the square norm of the residual
231 const Scalar result2 = residual.two_norm2();
232 using std::sqrt;
233 return sqrt(result2);
234 }
235
241 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
242 std::shared_ptr<SolutionVector> r)
243 {
244 jacobian_ = A;
245 residual_ = r;
246
247 setJacobianBuildMode(*jacobian_);
248 setJacobianPattern(*jacobian_);
249 setResidualSize(*residual_);
250 }
251
257 {
258 jacobian_ = std::make_shared<JacobianMatrix>();
259 residual_ = std::make_shared<SolutionVector>();
260
261 setJacobianBuildMode(*jacobian_);
262 setJacobianPattern(*jacobian_);
263 setResidualSize(*residual_);
264 }
265
270 {
271 using namespace Dune::Hybrid;
272 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
273 {
274 forEach(jac[i], [&](auto& jacBlock)
275 {
276 using BlockType = std::decay_t<decltype(jacBlock)>;
277 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
278 jacBlock.setBuildMode(BlockType::BuildMode::random);
279 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
280 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
281 });
282 });
283 }
284
289 {
290 using namespace Dune::Hybrid;
291 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
292 {
293 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
294 {
295 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
296 pattern.exportIdx(jac[domainI][domainJ]);
297 });
298 });
299 }
300
305 {
306 using namespace Dune::Hybrid;
307 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
308 { res[domainId].resize(this->numDofs(domainId)); });
309 }
310
315 {
316 using namespace Dune::Hybrid;
317 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
318 { this->gridVariables(domainId).update(curSol[domainId]); });
319 }
320
324 void resetTimeStep(const SolutionVector& curSol)
325 {
326 using namespace Dune::Hybrid;
327 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
328 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
329 }
330
332 template<std::size_t i>
333 std::size_t numDofs(Dune::index_constant<i> domainId) const
334 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
335
337 template<std::size_t i>
338 const auto& problem(Dune::index_constant<i> domainId) const
339 { return *std::get<domainId>(problemTuple_); }
340
342 template<std::size_t i>
343 const auto& gridGeometry(Dune::index_constant<i> domainId) const
344 { return *std::get<domainId>(gridGeometryTuple_); }
345
347 template<std::size_t i>
348 const auto& gridView(Dune::index_constant<i> domainId) const
349 { return gridGeometry(domainId).gridView(); }
350
352 template<std::size_t i>
353 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
354 { return *std::get<domainId>(gridVariablesTuple_); }
355
357 template<std::size_t i>
358 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
359 { return *std::get<domainId>(gridVariablesTuple_); }
360
363 { return *couplingManager_; }
364
367 { return *jacobian_; }
368
371 { return *residual_; }
372
374 const SolutionVector& prevSol() const
375 { return *prevSol_; }
376
381 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
382 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
383
389 { prevSol_ = &u; }
390
395 { return isStationaryProblem_; }
396
400 template<std::size_t i>
401 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
402 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
403
404protected:
406 std::shared_ptr<CouplingManager> couplingManager_;
407
408private:
409 // reset the residual vector to 0.0
410 void resetResidual_()
411 {
412 if(!residual_)
413 {
414 residual_ = std::make_shared<SolutionVector>();
415 setResidualSize(*residual_);
416 }
417
418 (*residual_) = 0.0;
419 }
420
421 // reset the jacobian vector to 0.0
422 void resetJacobian_()
423 {
424 if(!jacobian_)
425 {
426 jacobian_ = std::make_shared<JacobianMatrix>();
427 setJacobianBuildMode(*jacobian_);
428 setJacobianPattern(*jacobian_);
429 }
430
431 (*jacobian_) = 0.0;
432 }
433
434 // check if the assembler is in a correct state for assembly
435 void checkAssemblerState_() const
436 {
437 if (!isStationaryProblem_ && !prevSol_)
438 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
439
440 if (isStationaryProblem_ && prevSol_)
441 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
442 << " Did you forget to set the timeLoop to make this problem instationary?");
443 }
444
445 template<std::size_t i, class JacRow, class SubRes>
446 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
447 const SolutionVector& curSol)
448 {
449 assemble_(domainId, [&](const auto& element)
450 {
451 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
452 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
453 });
454 }
455
456 template<std::size_t i, class SubRes>
457 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
458 const SolutionVector& curSol)
459 {
460 assemble_(domainId, [&](const auto& element)
461 {
462 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
463 subDomainAssembler.assembleResidual(subRes);
464 });
465 }
466
473 template<std::size_t i, class AssembleElementFunc>
474 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
475 {
476 // let the local assembler add the element contributions
477 for (const auto& element : elements(gridView(domainId)))
478 assembleElement(element);
479 }
480
481 // get diagonal block pattern
482 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
483 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
484 Dune::index_constant<j> domainJ) const
485 {
486 const auto& gg = gridGeometry(domainI);
487 auto pattern = getJacobianPattern<isImplicit()>(gg);
488 couplingManager_->extendJacobianPattern(domainI, pattern);
489 return pattern;
490 }
491
492 // get coupling block pattern
493 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
494 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
495 Dune::index_constant<j> domainJ) const
496 {
497 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
498 domainI, gridGeometry(domainI),
499 domainJ, gridGeometry(domainJ));
500 }
501
503 ProblemTuple problemTuple_;
504
506 GridGeometryTuple gridGeometryTuple_;
507
509 GridVariablesTuple gridVariablesTuple_;
510
512 std::shared_ptr<const TimeLoop> timeLoop_;
513
515 const SolutionVector* prevSol_ = nullptr;
516
518 bool isStationaryProblem_;
519
521 std::shared_ptr<JacobianMatrix> jacobian_;
522 std::shared_ptr<SolutionVector> residual_;
523};
524
525} // end namespace Dumux
526
527#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.
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: common/properties.hh:77
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:62
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:333
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:314
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:85
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:86
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:82
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:366
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:93
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:207
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:84
CMType CouplingManager
Definition: multidomain/fvassembler.hh:88
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:200
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:324
Scalar residualNorm(const SolutionVector &curSol)
Definition: multidomain/fvassembler.hh:224
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:358
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:388
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:76
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:241
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:343
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:304
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:256
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:348
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:163
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:79
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:362
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:381
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:374
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:353
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:406
MDTraits Traits
Definition: multidomain/fvassembler.hh:67
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:269
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:73
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:143
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:401
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:338
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:184
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:394
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:69
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:288
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:370
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.
Manages the handling of time dependent problems.