version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
16
17#include <type_traits>
18#include <tuple>
19
20#include <dune/common/hybridutilities.hh>
21#include <dune/istl/matrixindexset.hh>
22
34
40#include "assemblerview.hh"
41
43
44namespace Dumux {
45
46namespace Grid::Capabilities {
47
48namespace Detail {
49// helper for multi-domain models
50template<class T, std::size_t... I>
51bool allGridsSupportsMultithreadingImpl(const T& gridGeometries, std::index_sequence<I...>)
52{
53 return (... && supportsMultithreading(std::get<I>(gridGeometries)->gridView()));
54}
55} // end namespace Detail
56
57// helper for multi-domain models (all grids have to support multithreading)
58template<class... GG>
59bool allGridsSupportsMultithreading(const std::tuple<GG...>& gridGeometries)
60{
61 return Detail::allGridsSupportsMultithreadingImpl<std::tuple<GG...>>(gridGeometries, std::make_index_sequence<sizeof...(GG)>());
62}
63
64} // end namespace Grid::Capabilities
65
71template<class CM>
72struct CouplingManagerSupportsMultithreadedAssembly : public std::false_type
73{};
74
84template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
86{
87 template<std::size_t id>
88 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
89
90public:
91 using Traits = MDTraits;
92
93 using Scalar = typename MDTraits::Scalar;
94
96 template<std::size_t id>
97 using LocalResidual = GetPropType<SubDomainTypeTag<id>, Properties::LocalResidual>;
98
99 template<std::size_t id>
100 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
101
102 template<std::size_t id>
103 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
104
105 template<std::size_t id>
106 using Problem = typename MDTraits::template SubDomain<id>::Problem;
107
108 using JacobianMatrix = typename MDTraits::JacobianMatrix;
109 using SolutionVector = typename MDTraits::SolutionVector;
110 using ResidualType = typename MDTraits::ResidualVector;
111
112 using CouplingManager = CMType;
113
117 static constexpr bool isImplicit()
118 { return useImplicitAssembly; }
119
120private:
121
122 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
123 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
124 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
125
127 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
128
129 template<std::size_t id>
130 using SubDomainAssemblerView = MultiDomainAssemblerSubDomainView<ThisType, id>;
131
132 template<class DiscretizationMethod, std::size_t id>
133 struct SubDomainAssemblerType;
134
135 template<std::size_t id>
136 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
137 {
138 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
139 };
140
141 template<std::size_t id>
142 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
143 {
144 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
145 };
146
147 template<std::size_t id, class DM>
148 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
149 {
150 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
151 };
152
153 template<std::size_t id>
154 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
155 {
156 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
157 };
158
159 template<std::size_t id>
160 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
161 {
162 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
163 };
164
165 template<std::size_t id>
166 using SubDomainAssembler = typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod, id>::type;
167
168public:
169
170
177 GridGeometryTuple gridGeometry,
178 GridVariablesTuple gridVariables,
179 std::shared_ptr<CouplingManager> couplingManager)
181 , problemTuple_(std::move(problem))
182 , gridGeometryTuple_(std::move(gridGeometry))
183 , gridVariablesTuple_(std::move(gridVariables))
184 , timeLoop_()
185 , isStationaryProblem_(true)
186 , warningIssued_(false)
187 {
188 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
189 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
190
194 && getParam<bool>("Assembly.Multithreading", true);
195
196 maybeComputeColors_();
197 }
198
205 GridGeometryTuple gridGeometry,
206 GridVariablesTuple gridVariables,
207 std::shared_ptr<CouplingManager> couplingManager,
208 std::shared_ptr<const TimeLoop> timeLoop,
209 const SolutionVector& prevSol)
211 , problemTuple_(std::move(problem))
212 , gridGeometryTuple_(std::move(gridGeometry))
213 , gridVariablesTuple_(std::move(gridVariables))
214 , timeLoop_(timeLoop)
215 , prevSol_(&prevSol)
216 , isStationaryProblem_(false)
217 , warningIssued_(false)
218 {
219 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
220
224 && getParam<bool>("Assembly.Multithreading", true);
225
226 maybeComputeColors_();
227 }
228
234 {
235 checkAssemblerState_();
236 resetJacobian_();
237 resetResidual_();
238
239 using namespace Dune::Hybrid;
240 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
241 {
242 auto& jacRow = (*jacobian_)[domainId];
243 auto& subRes = (*residual_)[domainId];
244 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
245
246 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
247 enforcePeriodicConstraints_(domainId, jacRow, subRes, *gridGeometry, curSol[domainId]);
248 });
249 }
250
253 {
254 resetResidual_();
255 assembleResidual(*residual_, curSol);
256 }
257
260 {
261 r = 0.0;
262
263 checkAssemblerState_();
264
265 // update the grid variables for the case of active caching
266 updateGridVariables(curSol);
267
268 using namespace Dune::Hybrid;
269 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
270 {
271 auto& subRes = r[domainId];
272 this->assembleResidual_(domainId, subRes, curSol);
273 });
274 }
275
281 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
282 std::shared_ptr<ResidualType> r)
283 {
284 jacobian_ = A;
285 residual_ = r;
286
287 setJacobianBuildMode(*jacobian_);
288 setJacobianPattern_(*jacobian_);
289 setResidualSize_(*residual_);
290 }
291
297 {
298 jacobian_ = std::make_shared<JacobianMatrix>();
299 residual_ = std::make_shared<ResidualType>();
300
301 setJacobianBuildMode(*jacobian_);
302 setJacobianPattern_(*jacobian_);
303 setResidualSize_(*residual_);
304 }
305
310 {
311 using namespace Dune::Hybrid;
312 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
313 {
314 forEach(jac[i], [&](auto& jacBlock)
315 {
316 using BlockType = std::decay_t<decltype(jacBlock)>;
317 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
318 jacBlock.setBuildMode(BlockType::BuildMode::random);
319 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
320 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
321 });
322 });
323 }
324
329 {
330 setJacobianPattern_(*jacobian_);
331 setResidualSize_(*residual_);
332 maybeComputeColors_();
333 }
334
339 {
340 using namespace Dune::Hybrid;
341 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
342 { this->gridVariables(domainId).update(curSol[domainId]); });
343 }
344
348 void resetTimeStep(const SolutionVector& curSol)
349 {
350 using namespace Dune::Hybrid;
351 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
352 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
353 }
354
356 template<std::size_t i>
357 std::size_t numDofs(Dune::index_constant<i> domainId) const
358 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
359
361 template<std::size_t i>
362 const auto& problem(Dune::index_constant<i> domainId) const
363 { return *std::get<domainId>(problemTuple_); }
364
366 template<std::size_t i>
367 const auto& gridGeometry(Dune::index_constant<i> domainId) const
368 { return *std::get<domainId>(gridGeometryTuple_); }
369
371 template<std::size_t i>
372 const auto& gridView(Dune::index_constant<i> domainId) const
373 { return gridGeometry(domainId).gridView(); }
374
376 template<std::size_t i>
377 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
378 { return *std::get<domainId>(gridVariablesTuple_); }
379
381 template<std::size_t i>
382 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
383 { return *std::get<domainId>(gridVariablesTuple_); }
384
387 { return *couplingManager_; }
388
391 { return *jacobian_; }
392
395 { return *residual_; }
396
398 const SolutionVector& prevSol() const
399 { return *prevSol_; }
400
405 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
406 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
407
413 { prevSol_ = &u; }
414
419 { return isStationaryProblem_; }
420
424 template<std::size_t i>
425 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
426 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
427
428protected:
430 std::shared_ptr<CouplingManager> couplingManager_;
431
432private:
436 void setJacobianPattern_(JacobianMatrix& jac) const
437 {
438 using namespace Dune::Hybrid;
439 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
440 {
441 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
442 {
443 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
444 pattern.exportIdx(jac[domainI][domainJ]);
445 });
446 });
447 }
448
452 void setResidualSize_(ResidualType& res) const
453 {
454 using namespace Dune::Hybrid;
455 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
456 { res[domainId].resize(this->numDofs(domainId)); });
457 }
458
459 // reset the residual vector to 0.0
460 void resetResidual_()
461 {
462 if(!residual_)
463 {
464 residual_ = std::make_shared<ResidualType>();
465 setResidualSize_(*residual_);
466 }
467
468 (*residual_) = 0.0;
469 }
470
471 // reset the jacobian vector to 0.0
472 void resetJacobian_()
473 {
474 if(!jacobian_)
475 {
476 jacobian_ = std::make_shared<JacobianMatrix>();
477 setJacobianBuildMode(*jacobian_);
478 setJacobianPattern_(*jacobian_);
479 }
480
481 (*jacobian_) = 0.0;
482 }
483
485 void maybeComputeColors_()
486 {
487 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
488 if (enableMultithreading_)
489 couplingManager_->computeColorsForAssembly();
490 }
491
492 // check if the assembler is in a correct state for assembly
493 void checkAssemblerState_() const
494 {
495 if (!isStationaryProblem_ && !prevSol_)
496 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
497
498 if (isStationaryProblem_ && prevSol_)
499 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
500 << " Did you forget to set the timeLoop to make this problem instationary?");
501 }
502
503 template<std::size_t i, class JacRow, class SubRes>
504 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
505 const SolutionVector& curSol)
506 {
507 assemble_(domainId, [&](const auto& element)
508 {
509 MultiDomainAssemblerSubDomainView view{*this, domainId};
510 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *couplingManager_);
511 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
512 });
513 }
514
515 template<std::size_t i, class SubRes>
516 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
517 const SolutionVector& curSol)
518 {
519 assemble_(domainId, [&](const auto& element)
520 {
521 MultiDomainAssemblerSubDomainView view{*this, domainId};
522 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *couplingManager_);
523 subDomainAssembler.assembleResidual(subRes);
524 });
525 }
526
532 template<std::size_t i, class AssembleElementFunc>
533 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
534 {
535 // a state that will be checked on all processes
536 bool succeeded = false;
537
538 // try assembling using the local assembly function
539 try
540 {
541 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
542 {
543 if (enableMultithreading_)
544 {
545 couplingManager_->assembleMultithreaded(
546 domainId, std::forward<AssembleElementFunc>(assembleElement)
547 );
548 return;
549 }
550 }
551
552 // fallback for coupling managers that don't support multithreaded assembly (yet)
553 // or if multithreaded assembly is disabled
554 // let the local assembler add the element contributions
555 for (const auto& element : elements(gridView(domainId)))
556 assembleElement(element);
557
558 // if we get here, everything worked well on this process
559 succeeded = true;
560 }
561 // throw exception if a problem occurred
562 catch (NumericalProblem &e)
563 {
564 std::cout << "rank " << gridView(domainId).comm().rank()
565 << " caught an exception while assembling:" << e.what()
566 << "\n";
567 succeeded = false;
568 }
569
570 // make sure everything worked well on all processes
571 if (gridView(domainId).comm().size() > 1)
572 succeeded = gridView(domainId).comm().min(succeeded);
573
574 // if not succeeded rethrow the error on all processes
575 if (!succeeded)
576 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
577 }
578
579 // get diagonal block pattern
580 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
581 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
582 Dune::index_constant<j> domainJ) const
583 {
584 const auto& gg = gridGeometry(domainI);
585 auto pattern = getJacobianPattern<isImplicit()>(gg);
586 couplingManager_->extendJacobianPattern(domainI, pattern);
587 return pattern;
588 }
589
590 // get coupling block pattern
591 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
592 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
593 Dune::index_constant<j> domainJ) const
594 {
595 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
596 domainI, gridGeometry(domainI),
597 domainJ, gridGeometry(domainJ));
598 }
599
600 // build periodic constraints into the system matrix
601 template<std::size_t i, class JacRow, class Res, class GG, class Sol>
602 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res, const GG& gridGeometry, const Sol& curSol)
603 {
604 if constexpr (Detail::hasPeriodicDofMap<GG>())
605 {
606 for (const auto& m : gridGeometry.periodicDofMap())
607 {
608 if (m.first < m.second)
609 {
610 auto& jac = jacRow[domainI];
611
612 // add the second row to the first
613 res[m.first] += res[m.second];
614
615 const auto end = jac[m.second].end();
616 for (auto it = jac[m.second].begin(); it != end; ++it)
617 jac[m.first][it.index()] += (*it);
618
619
620 // enforce the solution of the first periodic DOF to the second one
621 res[m.second] = curSol[m.second] - curSol[m.first];
622
623 // set derivatives accordingly in jacobian, i.e. id for m.second and -id for m.first
624 auto setMatrixBlock = [] (auto& matrixBlock, double diagValue)
625 {
626 for (int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
627 matrixBlock[eIdx][eIdx] = diagValue;
628 };
629
630 for (auto it = jac[m.second].begin(); it != end; ++it)
631 {
632 auto& matrixBlock = *it;
633 matrixBlock = 0.0;
634
635 assert(matrixBlock.N() == matrixBlock.M());
636 if(it.index() == m.second)
637 setMatrixBlock(matrixBlock, 1.0);
638
639 if(it.index() == m.first)
640 setMatrixBlock(matrixBlock, -1.0);
641
642 }
643
644 using namespace Dune::Hybrid;
645 forEach(makeIncompleteIntegerSequence<JacRow::size(), domainI>(), [&](const auto couplingDomainId)
646 {
647 auto& jacCoupling = jacRow[couplingDomainId];
648
649 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
650 jacCoupling[m.first][it.index()] += (*it);
651
652 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
653 (*it) = 0.0;
654 });
655 }
656 }
657 }
658 }
659
661 ProblemTuple problemTuple_;
662
664 GridGeometryTuple gridGeometryTuple_;
665
667 GridVariablesTuple gridVariablesTuple_;
668
670 std::shared_ptr<const TimeLoop> timeLoop_;
671
673 const SolutionVector* prevSol_ = nullptr;
674
676 bool isStationaryProblem_;
677
679 std::shared_ptr<JacobianMatrix> jacobian_;
680 std::shared_ptr<ResidualType> residual_;
681
683 mutable bool warningIssued_;
684
686 bool enableMultithreading_ = false;
687};
688
689} // end namespace Dumux
690
691#endif
Subdomain-specific views on multidomain assemblers.
Subdomain-specific view on a multidomain assembler. Allows retrieval of sub-domain specific objects w...
Definition: assemblerview.hh:31
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:86
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:357
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:338
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:328
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:109
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:106
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:390
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:117
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:259
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:108
CMType CouplingManager
Definition: multidomain/fvassembler.hh:112
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:252
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:348
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:382
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:412
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:100
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:367
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:296
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:372
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:103
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:386
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:405
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:398
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:377
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:430
MDTraits Traits
Definition: multidomain/fvassembler.hh:91
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:204
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:309
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:97
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:425
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:362
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:110
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:176
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:233
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:418
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:93
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:394
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< ResidualType > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: multidomain/fvassembler.hh:281
The cell-centered scheme multidomain local assembler.
Definition: multidomain/subdomaincclocalassembler.hh:270
The default time loop for instationary simulations.
Definition: common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Helper function to generate Jacobian pattern for multi domain models.
An enum class to define various differentiation methods available in order to compute the derivatives...
Some exceptions thrown in DuMux
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain proble...
Multithreading in Dumux.
bool allGridsSupportsMultithreadingImpl(const T &gridGeometries, std::index_sequence< I... >)
Definition: multistagemultidomainfvassembler.hh:52
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multistagemultidomainfvassembler.hh:60
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:74
Definition: adapt.hh:17
typename Detail::ConcatSeq< decltype(std::make_index_sequence< e >{}), e+1, decltype(std::make_index_sequence<(n > e) ?(n - e - 1) :0 >{})>::type makeIncompleteIntegerSequence
Definition: utility.hh:59
Provides a helper class for nonoverlapping decomposition.
Type traits to detect periodicity support.
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for...
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
Utilities for template meta programming.