version 3.11-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-FileCopyrightText: 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
39#include "assemblerview.hh"
40
42#if HAVE_DUMUX_OLD_STAGGERED
43#include <dumux/multidomain/subdomainstaggeredlocalassembler.hh>
44#endif
45
46namespace Dumux {
47
48template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM, bool implicit>
50
51namespace Grid::Capabilities {
52
53namespace Detail {
54// helper for multi-domain models
55template<class T, std::size_t... I>
56bool allGridsSupportsMultithreadingImpl(const T& gridGeometries, std::index_sequence<I...>)
57{
58 return (... && supportsMultithreading(std::get<I>(gridGeometries)->gridView()));
59}
60} // end namespace Detail
61
62// helper for multi-domain models (all grids have to support multithreading)
63template<class... GG>
64bool allGridsSupportsMultithreading(const std::tuple<GG...>& gridGeometries)
65{
66 return Detail::allGridsSupportsMultithreadingImpl<std::tuple<GG...>>(gridGeometries, std::make_index_sequence<sizeof...(GG)>());
67}
68
69} // end namespace Grid::Capabilities
70
71namespace Detail {
72
74template<class P>
75using SubProblemConstraintsDetector = decltype(std::declval<P>().constraints());
76
77template<class P>
78constexpr inline bool hasSubProblemGlobalConstraints()
79{ return Dune::Std::is_detected<SubProblemConstraintsDetector, P>::value; }
80
81} // end namespace Detail
82
88template<class CM>
89struct CouplingManagerSupportsMultithreadedAssembly : public std::false_type
90{};
91
101template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
103{
104 template<std::size_t id>
105 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
106
107public:
108 using Traits = MDTraits;
109
110 using Scalar = typename MDTraits::Scalar;
111
112 template<std::size_t id>
113 using LocalResidual = typename MDTraits::template SubDomain<id>::LocalResidual;
114
115 template<std::size_t id>
116 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
117
118 template<std::size_t id>
119 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
120
121 template<std::size_t id>
122 using Problem = typename MDTraits::template SubDomain<id>::Problem;
123
124 using JacobianMatrix = typename MDTraits::JacobianMatrix;
125 using SolutionVector = typename MDTraits::SolutionVector;
126 using ResidualType = typename MDTraits::ResidualVector;
127
128 using CouplingManager = CMType;
129
133 static constexpr bool isImplicit()
134 { return useImplicitAssembly; }
135
136private:
137
138 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
139 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
140 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
141
143 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
144
145 template<std::size_t id>
146 using SubDomainAssemblerView = MultiDomainAssemblerSubDomainView<ThisType, id>;
147
148 template<class DiscretizationMethod, std::size_t id>
149 struct SubDomainAssemblerType;
150
151 template<std::size_t id>
152 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
153 {
154 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
155 };
156
157 template<std::size_t id>
158 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
159 {
160 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
161 };
162
163 template<std::size_t id, class DM>
164 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
165 {
166 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
167 };
168
169 template<std::size_t id>
170 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
171 {
172 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
173 };
174
175 template<std::size_t id>
176 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
177 {
178 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod, isImplicit()>;
179 };
180
181 template<std::size_t id>
182 using SubDomainAssembler = typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod, id>::type;
183
184public:
185
186
193 GridGeometryTuple gridGeometry,
194 GridVariablesTuple gridVariables,
195 std::shared_ptr<CouplingManager> couplingManager)
197 , problemTuple_(std::move(problem))
198 , gridGeometryTuple_(std::move(gridGeometry))
199 , gridVariablesTuple_(std::move(gridVariables))
200 , timeLoop_()
201 , isStationaryProblem_(true)
202 , warningIssued_(false)
203 {
204 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
205 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
206
210 && getParam<bool>("Assembly.Multithreading", true);
211
212 maybeComputeColors_();
213 }
214
221 GridGeometryTuple gridGeometry,
222 GridVariablesTuple gridVariables,
223 std::shared_ptr<CouplingManager> couplingManager,
224 std::shared_ptr<const TimeLoop> timeLoop,
225 const SolutionVector& prevSol)
227 , problemTuple_(std::move(problem))
228 , gridGeometryTuple_(std::move(gridGeometry))
229 , gridVariablesTuple_(std::move(gridVariables))
230 , timeLoop_(timeLoop)
231 , prevSol_(&prevSol)
232 , isStationaryProblem_(false)
233 , warningIssued_(false)
234 {
235 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
236
240 && getParam<bool>("Assembly.Multithreading", true);
241
242 maybeComputeColors_();
243 }
244
250 {
251 checkAssemblerState_();
252 resetJacobian_();
253 resetResidual_();
254
255 using namespace Dune::Hybrid;
256 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
257 {
258 auto& jacRow = (*jacobian_)[domainId];
259 auto& subRes = (*residual_)[domainId];
260 const auto& subSol = curSol[domainId];
261 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
262
263 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
264 enforcePeriodicConstraints_(domainId, jacRow, subRes, *gridGeometry, subSol);
265 enforceProblemConstraints_(domainId, jacRow, subRes, *gridGeometry, subSol);
266 });
267 }
268
271 {
272 resetResidual_();
273 assembleResidual(*residual_, curSol);
274 }
275
278 {
279 r = 0.0;
280
281 checkAssemblerState_();
282
283 // update the grid variables for the case of active caching
284 updateGridVariables(curSol);
285
286 using namespace Dune::Hybrid;
287 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
288 {
289 auto& subRes = r[domainId];
290 const auto& subSol = curSol[domainId];
291 this->assembleResidual_(domainId, subRes, curSol);
292
293 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
294 enforcePeriodicConstraints_(domainId, subRes, *gridGeometry, subSol);
295 enforceProblemConstraints_(domainId, subRes, *gridGeometry, subSol);
296 });
297 }
298
304 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
305 std::shared_ptr<ResidualType> r)
306 {
307 jacobian_ = A;
308 residual_ = r;
309
310 setJacobianBuildMode(*jacobian_);
311 setJacobianPattern_(*jacobian_);
312 setResidualSize_(*residual_);
313 }
314
320 {
321 jacobian_ = std::make_shared<JacobianMatrix>();
322 residual_ = std::make_shared<ResidualType>();
323
324 setJacobianBuildMode(*jacobian_);
325 setJacobianPattern_(*jacobian_);
326 setResidualSize_(*residual_);
327 }
328
333 {
334 using namespace Dune::Hybrid;
335 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
336 {
337 forEach(jac[i], [&](auto& jacBlock)
338 {
339 using BlockType = std::decay_t<decltype(jacBlock)>;
340 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
341 jacBlock.setBuildMode(BlockType::BuildMode::random);
342 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
343 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
344 });
345 });
346 }
347
352 {
353 setJacobianPattern_(*jacobian_);
354 setResidualSize_(*residual_);
355 maybeComputeColors_();
356 }
357
362 {
363 using namespace Dune::Hybrid;
364 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
365 { this->gridVariables(domainId).update(curSol[domainId]); });
366 }
367
371 void resetTimeStep(const SolutionVector& curSol)
372 {
373 using namespace Dune::Hybrid;
374 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
375 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
376 }
377
379 template<std::size_t i>
380 std::size_t numDofs(Dune::index_constant<i> domainId) const
381 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
382
384 template<std::size_t i>
385 const auto& problem(Dune::index_constant<i> domainId) const
386 { return *std::get<domainId>(problemTuple_); }
387
389 template<std::size_t i>
390 const auto& gridGeometry(Dune::index_constant<i> domainId) const
391 { return *std::get<domainId>(gridGeometryTuple_); }
392
394 template<std::size_t i>
395 const auto& gridView(Dune::index_constant<i> domainId) const
396 { return gridGeometry(domainId).gridView(); }
397
399 template<std::size_t i>
400 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
401 { return *std::get<domainId>(gridVariablesTuple_); }
402
404 template<std::size_t i>
405 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
406 { return *std::get<domainId>(gridVariablesTuple_); }
407
410 { return *couplingManager_; }
411
414 { return *jacobian_; }
415
418 { return *residual_; }
419
421 const SolutionVector& prevSol() const
422 { return *prevSol_; }
423
428 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
429 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
430
436 { prevSol_ = &u; }
437
442 { return isStationaryProblem_; }
443
447 template<std::size_t i>
448 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
449 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
450
451protected:
453 std::shared_ptr<CouplingManager> couplingManager_;
454
455private:
459 void setJacobianPattern_(JacobianMatrix& jac) const
460 {
461 using namespace Dune::Hybrid;
462 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
463 {
464 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
465 {
466 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
467 pattern.exportIdx(jac[domainI][domainJ]);
468 });
469 });
470 }
471
475 void setResidualSize_(ResidualType& res) const
476 {
477 using namespace Dune::Hybrid;
478 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
479 { res[domainId].resize(this->numDofs(domainId)); });
480 }
481
482 // reset the residual vector to 0.0
483 void resetResidual_()
484 {
485 if(!residual_)
486 {
487 residual_ = std::make_shared<ResidualType>();
488 setResidualSize_(*residual_);
489 }
490
491 (*residual_) = 0.0;
492 }
493
494 // reset the jacobian vector to 0.0
495 void resetJacobian_()
496 {
497 if(!jacobian_)
498 {
499 jacobian_ = std::make_shared<JacobianMatrix>();
500 setJacobianBuildMode(*jacobian_);
501 setJacobianPattern_(*jacobian_);
502 }
503
504 (*jacobian_) = 0.0;
505 }
506
508 void maybeComputeColors_()
509 {
510 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
511 if (enableMultithreading_)
512 couplingManager_->computeColorsForAssembly();
513 }
514
515 // check if the assembler is in a correct state for assembly
516 void checkAssemblerState_() const
517 {
518 if (!isStationaryProblem_ && !prevSol_)
519 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
520
521 if (isStationaryProblem_ && prevSol_)
522 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
523 << " Did you forget to set the timeLoop to make this problem instationary?");
524 }
525
526 template<std::size_t i, class JacRow, class SubRes>
527 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
528 const SolutionVector& curSol)
529 {
530 assemble_(domainId, [&](const auto& element)
531 {
532 MultiDomainAssemblerSubDomainView view{*this, domainId};
533 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *couplingManager_);
534 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
535 });
536 }
537
538 template<std::size_t i, class SubRes>
539 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
540 const SolutionVector& curSol)
541 {
542 assemble_(domainId, [&](const auto& element)
543 {
544 MultiDomainAssemblerSubDomainView view{*this, domainId};
545 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *couplingManager_);
546 subDomainAssembler.assembleResidual(subRes);
547 });
548 }
549
555 template<std::size_t i, class AssembleElementFunc>
556 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
557 {
558 // a state that will be checked on all processes
559 bool succeeded = false;
560
561 // try assembling using the local assembly function
562 try
563 {
564 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
565 {
566 if (enableMultithreading_)
567 {
568 couplingManager_->assembleMultithreaded(
569 domainId, std::forward<AssembleElementFunc>(assembleElement)
570 );
571 return;
572 }
573 }
574
575 // fallback for coupling managers that don't support multithreaded assembly (yet)
576 // or if multithreaded assembly is disabled
577 // let the local assembler add the element contributions
578 for (const auto& element : elements(gridView(domainId)))
579 assembleElement(element);
580
581 // if we get here, everything worked well on this process
582 succeeded = true;
583 }
584 // throw exception if a problem occurred
585 catch (NumericalProblem &e)
586 {
587 std::cout << "rank " << gridView(domainId).comm().rank()
588 << " caught an exception while assembling:" << e.what()
589 << "\n";
590 succeeded = false;
591 }
592
593 // make sure everything worked well on all processes
594 if (gridView(domainId).comm().size() > 1)
595 succeeded = gridView(domainId).comm().min(succeeded);
596
597 // if not succeeded rethrow the error on all processes
598 if (!succeeded)
599 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
600 }
601
602 // get diagonal block pattern
603 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
604 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
605 Dune::index_constant<j> domainJ) const
606 {
607 const auto& gg = gridGeometry(domainI);
608 auto pattern = getJacobianPattern<isImplicit()>(gg);
609 couplingManager_->extendJacobianPattern(domainI, pattern);
610 return pattern;
611 }
612
613 // get coupling block pattern
614 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
615 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
616 Dune::index_constant<j> domainJ) const
617 {
618 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
619 domainI, gridGeometry(domainI),
620 domainJ, gridGeometry(domainJ));
621 }
622
623 // build periodic constraints into the system matrix
624 template<std::size_t i, class JacRow, class Res, class GG, class Sol>
625 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res, const GG& gridGeometry, const Sol& curSol) const
626 {
627 if constexpr (Detail::hasPeriodicDofMap<GG>())
628 {
629 for (const auto& m : gridGeometry.periodicDofMap())
630 {
631 if (m.first < m.second)
632 {
633 auto& jac = jacRow[domainI];
634
635 // add the second row to the first
636 res[m.first] += res[m.second];
637
638 const auto end = jac[m.second].end();
639 for (auto it = jac[m.second].begin(); it != end; ++it)
640 jac[m.first][it.index()] += (*it);
641
642
643 // enforce the solution of the first periodic DOF to the second one
644 res[m.second] = curSol[m.second] - curSol[m.first];
645
646 // set derivatives accordingly in jacobian, i.e. id for m.second and -id for m.first
647 auto setMatrixBlock = [] (auto& matrixBlock, double diagValue)
648 {
649 for (int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
650 matrixBlock[eIdx][eIdx] = diagValue;
651 };
652
653 for (auto it = jac[m.second].begin(); it != end; ++it)
654 {
655 auto& matrixBlock = *it;
656 matrixBlock = 0.0;
657
658 assert(matrixBlock.N() == matrixBlock.M());
659 if(it.index() == m.second)
660 setMatrixBlock(matrixBlock, 1.0);
661
662 if(it.index() == m.first)
663 setMatrixBlock(matrixBlock, -1.0);
664
665 }
666
667 using namespace Dune::Hybrid;
668 forEach(makeIncompleteIntegerSequence<JacRow::size(), domainI>(), [&](const auto couplingDomainId)
669 {
670 auto& jacCoupling = jacRow[couplingDomainId];
671
672 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
673 jacCoupling[m.first][it.index()] += (*it);
674
675 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
676 (*it) = 0.0;
677 });
678 }
679 }
680 }
681 }
682
683 // enforce global constraints into the system matrix
684 template<std::size_t i, class JacRow, class Res, class GG, class Sol>
685 void enforceProblemConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res, const GG& gridGeometry, const Sol& curSol) const
686 {
687 if constexpr (Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
688 {
689 auto& jac = jacRow[domainI];
690 auto applyDirichletConstraint = [&] (const auto& dofIdx,
691 const auto& values,
692 const auto eqIdx,
693 const auto pvIdx)
694 {
695 (res)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
696
697 auto& row = jac[dofIdx];
698 for (auto col = row.begin(); col != row.end(); ++col)
699 row[col.index()][eqIdx] = 0.0;
700
701 jac[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
702
703 using namespace Dune::Hybrid;
704 forEach(makeIncompleteIntegerSequence<JacRow::size(), domainI>(), [&](const auto couplingDomainId)
705 {
706 auto& jacCoupling = jacRow[couplingDomainId];
707
708 auto& rowCoupling = jacCoupling[dofIdx];
709 for (auto colCoupling = rowCoupling.begin(); colCoupling != rowCoupling.end(); ++colCoupling)
710 rowCoupling[colCoupling.index()][eqIdx] = 0.0;
711 });
712
713 };
714
715 for (const auto& constraintData : this->problem(domainI).constraints())
716 {
717 const auto& constraintInfo = constraintData.constraintInfo();
718 const auto& values = constraintData.values();
719 const auto dofIdx = constraintData.dofIndex();
720 // set the constraint in residual and jacobian
721 for (int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
722 {
723 if (constraintInfo.isConstraintEquation(eqIdx))
724 {
725 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
726 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
727 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
728 }
729 }
730 }
731 }
732 }
733
734 // build periodic constraints in residual only
735 template<std::size_t i, class Res, class GG, class Sol>
736 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, Res& res, const GG& gridGeometry, const Sol& curSol) const
737 {
738 if constexpr (Detail::hasPeriodicDofMap<GG>())
739 {
740 for (const auto& m : gridGeometry.periodicDofMap())
741 {
742 if (m.first < m.second)
743 {
744 // add the second row to the first
745 res[m.first] += res[m.second];
746
747 // enforce the solution of the first periodic DOF to the second one
748 res[m.second] = curSol[m.second] - curSol[m.first];
749 }
750 }
751 }
752 }
753
754 // enforce global constraints in residual only
755 template<std::size_t i, class Res, class GG, class Sol>
756 void enforceProblemConstraints_(Dune::index_constant<i> domainI, Res& res, const GG& gridGeometry, const Sol& curSol) const
757 {
758 if constexpr (Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
759 {
760 for (const auto& constraintData : this->problem(domainI).constraints())
761 {
762 const auto& constraintInfo = constraintData.constraintInfo();
763 const auto& values = constraintData.values();
764 const auto dofIdx = constraintData.dofIndex();
765 for (int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
766 {
767 if (constraintInfo.isConstraintEquation(eqIdx))
768 {
769 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
770 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
771 res[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
772 }
773 }
774 }
775 }
776 }
777
779 ProblemTuple problemTuple_;
780
782 GridGeometryTuple gridGeometryTuple_;
783
785 GridVariablesTuple gridVariablesTuple_;
786
788 std::shared_ptr<const TimeLoop> timeLoop_;
789
791 const SolutionVector* prevSol_ = nullptr;
792
794 bool isStationaryProblem_;
795
797 std::shared_ptr<JacobianMatrix> jacobian_;
798 std::shared_ptr<ResidualType> residual_;
799
801 mutable bool warningIssued_;
802
804 bool enableMultithreading_ = false;
805};
806
807} // end namespace Dumux
808
809#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:103
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:380
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:361
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:351
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:125
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:122
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:413
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:133
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:277
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:124
CMType CouplingManager
Definition: multidomain/fvassembler.hh:128
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:270
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:371
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:405
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:435
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:116
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:390
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:319
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:395
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:119
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:409
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:428
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:421
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:400
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:453
MDTraits Traits
Definition: multidomain/fvassembler.hh:108
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:220
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:332
typename MDTraits::template SubDomain< id >::LocalResidual LocalResidual
Definition: multidomain/fvassembler.hh:113
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:448
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:385
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:126
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:192
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:249
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:441
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:110
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:417
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:304
The cell-centered scheme multidomain local assembler.
Definition: multidomain/subdomaincclocalassembler.hh:270
Definition: multidomain/fvassembler.hh:49
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
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.
decltype(std::declval< P >().constraints()) SubProblemConstraintsDetector
helper struct detecting if sub-problem has a constraints() function
Definition: multidomain/fvassembler.hh:75
constexpr bool hasSubProblemGlobalConstraints()
Definition: multidomain/fvassembler.hh:78
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:34
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...
Utilities for template meta programming.