3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Loading...
Searching...
No Matches
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
40
45
46namespace Dumux {
47
57template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
59{
60 template<std::size_t id>
61 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
62
63public:
64 using Traits = MDTraits;
65
66 using Scalar = typename MDTraits::Scalar;
67
69 template<std::size_t id>
71
72 template<std::size_t id>
73 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
74
75 template<std::size_t id>
76 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
77
78 template<std::size_t id>
79 using FVGridGeometry [[deprecated("Use GridGeometry instead. FVGridGeometry will be removed after 3.1!")]] = GridGeometry<id>;
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
162 [[deprecated("Please use constructor taking the previous solution instead. Will be removed after release 3.2!")]]
164 GridGeometryTuple&& gridGeometry,
165 GridVariablesTuple&& gridVariables,
166 std::shared_ptr<CouplingManager> couplingManager,
167 std::shared_ptr<const TimeLoop> timeLoop)
169 , problemTuple_(problem)
170 , gridGeometryTuple_(gridGeometry)
171 , gridVariablesTuple_(gridVariables)
172 , timeLoop_(timeLoop)
173 , isStationaryProblem_(false)
174 {
175 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
176 }
177
184 GridGeometryTuple&& gridGeometry,
185 GridVariablesTuple&& gridVariables,
186 std::shared_ptr<CouplingManager> couplingManager,
187 std::shared_ptr<const TimeLoop> timeLoop,
188 const SolutionVector& prevSol)
190 , problemTuple_(problem)
191 , gridGeometryTuple_(gridGeometry)
192 , gridVariablesTuple_(gridVariables)
193 , timeLoop_(timeLoop)
194 , prevSol_(&prevSol)
195 , isStationaryProblem_(false)
196 {
197 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
198 }
199
205 {
206 checkAssemblerState_();
207 resetJacobian_();
208 resetResidual_();
209
210 using namespace Dune::Hybrid;
211 forEach(integralRange(Dune::Hybrid::size(*jacobian_)), [&](const auto domainId)
212 {
213 auto& jacRow = (*jacobian_)[domainId];
214 auto& subRes = (*residual_)[domainId];
215 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
216 });
217 }
218
221 {
222 resetResidual_();
223 assembleResidual(*residual_, curSol);
224 }
225
228 {
229 checkAssemblerState_();
230
231 // update the grid variables for the case of active caching
232 updateGridVariables(curSol);
233
234 using namespace Dune::Hybrid;
235 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
236 {
237 auto& subRes = r[domainId];
238 this->assembleResidual_(domainId, subRes, curSol);
239 });
240 }
241
245 {
248 assembleResidual(residual, curSol);
249
250 // calculate the square norm of the residual
251 const Scalar result2 = residual.two_norm2();
252 using std::sqrt;
253 return sqrt(result2);
254 }
255
261 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
262 std::shared_ptr<SolutionVector> r)
263 {
264 jacobian_ = A;
265 residual_ = r;
266
267 setJacobianBuildMode(*jacobian_);
268 setJacobianPattern(*jacobian_);
269 setResidualSize(*residual_);
270 }
271
277 {
278 jacobian_ = std::make_shared<JacobianMatrix>();
279 residual_ = std::make_shared<SolutionVector>();
280
281 setJacobianBuildMode(*jacobian_);
282 setJacobianPattern(*jacobian_);
283 setResidualSize(*residual_);
284 }
285
290 {
291 using namespace Dune::Hybrid;
292 forEach(jac, [](auto& jacRow)
293 {
294 forEach(jacRow, [](auto& jacBlock)
295 {
296 using BlockType = std::decay_t<decltype(jacBlock)>;
297 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
298 jacBlock.setBuildMode(BlockType::BuildMode::random);
299 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
300 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
301 });
302 });
303 }
304
309 {
310 using namespace Dune::Hybrid;
311 forEach(integralRange(Dune::Hybrid::size(jac)), [&](const auto domainI)
312 {
313 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
314 {
315 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
316 pattern.exportIdx(jac[domainI][domainJ]);
317 });
318 });
319 }
320
325 {
326 using namespace Dune::Hybrid;
327 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
328 { res[domainId].resize(this->numDofs(domainId)); });
329 }
330
335 {
336 using namespace Dune::Hybrid;
337 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
338 { this->gridVariables(domainId).update(curSol[domainId]); });
339 }
340
344 void resetTimeStep(const SolutionVector& curSol)
345 {
346 using namespace Dune::Hybrid;
347 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
348 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
349 }
350
352 template<std::size_t i>
353 std::size_t numDofs(Dune::index_constant<i> domainId) const
354 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
355
357 template<std::size_t i>
358 const auto& problem(Dune::index_constant<i> domainId) const
359 { return *std::get<domainId>(problemTuple_); }
360
362 template<std::size_t i>
363 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
364 const auto& fvGridGeometry(Dune::index_constant<i> domainId) const
365 { return gridGeometry(domainId); }
366
368 template<std::size_t i>
369 const auto& gridGeometry(Dune::index_constant<i> domainId) const
370 { return *std::get<domainId>(gridGeometryTuple_); }
371
373 template<std::size_t i>
374 const auto& gridView(Dune::index_constant<i> domainId) const
375 { return gridGeometry(domainId).gridView(); }
376
378 template<std::size_t i>
379 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
380 { return *std::get<domainId>(gridVariablesTuple_); }
381
383 template<std::size_t i>
384 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
385 { return *std::get<domainId>(gridVariablesTuple_); }
386
389 { return *couplingManager_; }
390
393 { return *jacobian_; }
394
397 { return *residual_; }
398
400 const SolutionVector& prevSol() const
401 { return *prevSol_; }
402
407 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
408 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
409
415 { prevSol_ = &u; }
416
421 { return isStationaryProblem_; }
422
426 template<std::size_t i>
427 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
428 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
429
430protected:
432 std::shared_ptr<CouplingManager> couplingManager_;
433
434private:
435 // reset the residual vector to 0.0
436 void resetResidual_()
437 {
438 if(!residual_)
439 {
440 residual_ = std::make_shared<SolutionVector>();
441 setResidualSize(*residual_);
442 }
443
444 (*residual_) = 0.0;
445 }
446
447 // reset the jacobian vector to 0.0
448 void resetJacobian_()
449 {
450 if(!jacobian_)
451 {
452 jacobian_ = std::make_shared<JacobianMatrix>();
453 setJacobianBuildMode(*jacobian_);
454 setJacobianPattern(*jacobian_);
455 }
456
457 (*jacobian_) = 0.0;
458 }
459
460 // check if the assembler is in a correct state for assembly
461 void checkAssemblerState_() const
462 {
463 if (!isStationaryProblem_ && !prevSol_)
464 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
465
466 if (isStationaryProblem_ && prevSol_)
467 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
468 << " Did you forget to set the timeLoop to make this problem instationary?");
469 }
470
471 template<std::size_t i, class JacRow, class SubRes>
472 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
473 const SolutionVector& curSol)
474 {
475 assemble_(domainId, [&](const auto& element)
476 {
477 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
478 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
479 });
480 }
481
482 template<std::size_t i, class SubRes>
483 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
484 const SolutionVector& curSol)
485 {
486 assemble_(domainId, [&](const auto& element)
487 {
488 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
489 subDomainAssembler.assembleResidual(subRes);
490 });
491 }
492
499 template<std::size_t i, class AssembleElementFunc>
500 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
501 {
502 // let the local assembler add the element contributions
503 for (const auto& element : elements(gridView(domainId)))
504 assembleElement(element);
505 }
506
507 // get diagonal block pattern
508 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
509 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
510 Dune::index_constant<j> domainJ) const
511 {
512 const auto& gg = gridGeometry(domainI);
513 auto pattern = getJacobianPattern<isImplicit()>(gg);
514 couplingManager_->extendJacobianPattern(domainI, pattern);
515 return pattern;
516 }
517
518 // get coupling block pattern
519 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
520 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
521 Dune::index_constant<j> domainJ) const
522 {
524 domainI, gridGeometry(domainI),
525 domainJ, gridGeometry(domainJ));
526 }
527
529 ProblemTuple problemTuple_;
530
532 GridGeometryTuple gridGeometryTuple_;
533
535 GridVariablesTuple gridVariablesTuple_;
536
538 std::shared_ptr<const TimeLoop> timeLoop_;
539
541 const SolutionVector* prevSol_ = nullptr;
542
544 bool isStationaryProblem_;
545
547 std::shared_ptr<JacobianMatrix> jacobian_;
548 std::shared_ptr<SolutionVector> residual_;
549};
550
551} // end namespace Dumux
552
553#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.
Manages the handling of time dependent problems.
The available discretization methods in Dumux.
Helper function to generate Jacobian pattern for multi domain models.
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...
A multidomain assembler for Jacobian and residual contribution per element (staggered method).
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ box
Definition method.hh:38
@ ccmpfa
Definition method.hh:38
@ staggered
Definition method.hh:38
@ cctpfa
Definition method.hh:38
Dune::MatrixIndexSet getJacobianPattern(const GridGeometry &gridGeometry)
Helper function to generate Jacobian pattern for the box method.
Definition jacobianpattern.hh:39
Dune::MatrixIndexSet getCouplingJacobianPattern(const CouplingManager &couplingManager, Dune::index_constant< i > domainI, const GridGeometryI &gridGeometryI, Dune::index_constant< j > domainJ, const GridGeometryJ &gridGeometryJ)
Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for cell-centered schemes...
Definition couplingjacobianpattern.hh:42
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
const Scalar PengRobinsonMixture< Scalar, StaticParameters >::u
Definition pengrobinsonmixture.hh:167
Definition common/properties.hh:91
Manages the handling of time dependent problems.
Definition timeloop.hh:65
The default time loop for instationary simulations.
Definition timeloop.hh:110
Definition multidomain/couplingmanager.hh:46
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition multidomain/fvassembler.hh:59
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition multidomain/fvassembler.hh:353
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition multidomain/fvassembler.hh:334
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::SolutionVector SolutionVector
Definition multidomain/fvassembler.hh:85
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::Problem Problem
Definition multidomain/fvassembler.hh:82
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition multidomain/fvassembler.hh:392
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:227
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::JacobianMatrix JacobianMatrix
Definition multidomain/fvassembler.hh:84
StaggeredCouplingManager< StaggeredMultiDomainTraits< TypeTag, TypeTag > > CouplingManager
Definition multidomain/fvassembler.hh:88
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition multidomain/fvassembler.hh:220
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition multidomain/fvassembler.hh:344
Scalar residualNorm(const SolutionVector &curSol)
Definition multidomain/fvassembler.hh:244
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition multidomain/fvassembler.hh:384
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:414
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::GridVariables GridVariables
Definition multidomain/fvassembler.hh:73
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:261
const auto & fvGridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition multidomain/fvassembler.hh:364
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition multidomain/fvassembler.hh:324
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition multidomain/fvassembler.hh:276
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition multidomain/fvassembler.hh:374
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:183
typename StaggeredMultiDomainTraits< TypeTag, TypeTag >::template SubDomain< id >::GridGeometry GridGeometry
Definition multidomain/fvassembler.hh:76
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop)
The constructor for instationary problems.
Definition multidomain/fvassembler.hh:163
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition multidomain/fvassembler.hh:407
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition multidomain/fvassembler.hh:289
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:427
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition multidomain/fvassembler.hh:204
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition multidomain/fvassembler.hh:420
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition multidomain/fvassembler.hh:308
The cell-centered scheme multidomain local assembler.
Definition subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.