3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/fvassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
28
29#include <type_traits>
30#include <tuple>
31
32#include <dune/common/hybridutilities.hh>
33#include <dune/istl/matrixindexset.hh>
34
42
48
50
51namespace Dumux {
52
62template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
64{
65 template<std::size_t id>
66 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
67
68public:
69 using Traits = MDTraits;
70
71 using Scalar = typename MDTraits::Scalar;
72
74 template<std::size_t id>
76
77 template<std::size_t id>
78 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
79
80 template<std::size_t id>
81 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
82
83 template<std::size_t id>
84 using Problem = typename MDTraits::template SubDomain<id>::Problem;
85
86 using JacobianMatrix = typename MDTraits::JacobianMatrix;
87 using SolutionVector = typename MDTraits::SolutionVector;
89
90 using CouplingManager = CMType;
91
95 static constexpr bool isImplicit()
96 { return useImplicitAssembly; }
97
98private:
99
100 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
101 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
102 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
103
105 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
106
107 template<class DiscretizationMethod, std::size_t id>
108 struct SubDomainAssemblerType;
109
110 template<std::size_t id>
111 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
112 {
113 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
114 };
115
116 template<std::size_t id>
117 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
118 {
119 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
120 };
121
122 template<std::size_t id>
123 struct SubDomainAssemblerType<DiscretizationMethods::Box, id>
124 {
125 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
126 };
127
128 template<std::size_t id>
129 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
130 {
131 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
132 };
133
134 template<std::size_t id>
135 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
136 {
137 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
138 };
139
140 template<std::size_t id>
141 using SubDomainAssembler = typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod, id>::type;
142
143public:
144
145
152 GridGeometryTuple gridGeometry,
153 GridVariablesTuple gridVariables,
154 std::shared_ptr<CouplingManager> couplingManager)
156 , problemTuple_(std::move(problem))
157 , gridGeometryTuple_(std::move(gridGeometry))
158 , gridVariablesTuple_(std::move(gridVariables))
159 , timeLoop_()
160 , isStationaryProblem_(true)
161 , warningIssued_(false)
162 {
163 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
164 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
165 }
166
173 GridGeometryTuple gridGeometry,
174 GridVariablesTuple gridVariables,
175 std::shared_ptr<CouplingManager> couplingManager,
176 std::shared_ptr<const TimeLoop> timeLoop,
177 const SolutionVector& prevSol)
179 , problemTuple_(std::move(problem))
180 , gridGeometryTuple_(std::move(gridGeometry))
181 , gridVariablesTuple_(std::move(gridVariables))
182 , timeLoop_(timeLoop)
183 , prevSol_(&prevSol)
184 , isStationaryProblem_(false)
185 , warningIssued_(false)
186 {
187 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
188 }
189
195 {
196 checkAssemblerState_();
197 resetJacobian_();
198 resetResidual_();
199
200 using namespace Dune::Hybrid;
201 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
202 {
203 auto& jacRow = (*jacobian_)[domainId];
204 auto& subRes = (*residual_)[domainId];
205 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
206
207 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
208 enforcePeriodicConstraints_(domainId, jacRow, subRes, *gridGeometry, curSol[domainId]);
209 });
210 }
211
214 {
215 resetResidual_();
216 assembleResidual(*residual_, curSol);
217 }
218
221 {
222 r = 0.0;
223
224 checkAssemblerState_();
225
226 // update the grid variables for the case of active caching
227 updateGridVariables(curSol);
228
229 using namespace Dune::Hybrid;
230 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
231 {
232 auto& subRes = r[domainId];
233 this->assembleResidual_(domainId, subRes, curSol);
234 });
235 }
236
239 {
242 assembleResidual(residual, curSol);
243
244 // calculate the squared norm of the residual
245 Scalar resultSquared = 0.0;
246
247 // for box communicate the residual with the neighboring processes
248 using namespace Dune::Hybrid;
249 forEach(integralRange(Dune::Hybrid::size(residual)), [&](const auto domainId)
250 {
251 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
252 const auto& gridView = gridGeometry->gridView();
253
254 if (gridView.comm().size() > 1 && gridView.overlapSize(0) == 0)
255 {
257 {
258 using GV = typename GridGeometry<domainId>::GridView;
259 using DM = typename GridGeometry<domainId>::VertexMapper;
261
262 PVHelper vectorHelper(gridView, gridGeometry->vertexMapper());
263
264 vectorHelper.makeNonOverlappingConsistent(residual[domainId]);
265 }
266 }
267 else if (!warningIssued_)
268 {
269 if (gridView.comm().size() > 1 && gridView.comm().rank() == 0)
270 std::cout << "\nWarning: norm calculation adds entries corresponding to\n"
271 << "overlapping entities multiple times. Please use the norm\n"
272 << "function provided by a linear solver instead." << std::endl;
273
274 warningIssued_ = true;
275 }
276
277 Scalar localNormSquared = residual[domainId].two_norm2();
278
279 if (gridView.comm().size() > 1)
280 {
281 localNormSquared = gridView.comm().sum(localNormSquared);
282 }
283
284 resultSquared += localNormSquared;
285 });
286
287 using std::sqrt;
288 return sqrt(resultSquared);
289 }
290
296 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
297 std::shared_ptr<SolutionVector> r)
298 {
299 jacobian_ = A;
300 residual_ = r;
301
302 setJacobianBuildMode(*jacobian_);
303 setJacobianPattern(*jacobian_);
304 setResidualSize(*residual_);
305 }
306
312 {
313 jacobian_ = std::make_shared<JacobianMatrix>();
314 residual_ = std::make_shared<SolutionVector>();
315
316 setJacobianBuildMode(*jacobian_);
317 setJacobianPattern(*jacobian_);
318 setResidualSize(*residual_);
319 }
320
325 {
326 using namespace Dune::Hybrid;
327 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
328 {
329 forEach(jac[i], [&](auto& jacBlock)
330 {
331 using BlockType = std::decay_t<decltype(jacBlock)>;
332 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
333 jacBlock.setBuildMode(BlockType::BuildMode::random);
334 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
335 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
336 });
337 });
338 }
339
344 {
345 using namespace Dune::Hybrid;
346 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
347 {
348 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
349 {
350 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
351 pattern.exportIdx(jac[domainI][domainJ]);
352 });
353 });
354 }
355
360 {
361 using namespace Dune::Hybrid;
362 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
363 { res[domainId].resize(this->numDofs(domainId)); });
364 }
365
370 {
371 using namespace Dune::Hybrid;
372 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
373 { this->gridVariables(domainId).update(curSol[domainId]); });
374 }
375
379 void resetTimeStep(const SolutionVector& curSol)
380 {
381 using namespace Dune::Hybrid;
382 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
383 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
384 }
385
387 template<std::size_t i>
388 std::size_t numDofs(Dune::index_constant<i> domainId) const
389 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
390
392 template<std::size_t i>
393 const auto& problem(Dune::index_constant<i> domainId) const
394 { return *std::get<domainId>(problemTuple_); }
395
397 template<std::size_t i>
398 const auto& gridGeometry(Dune::index_constant<i> domainId) const
399 { return *std::get<domainId>(gridGeometryTuple_); }
400
402 template<std::size_t i>
403 const auto& gridView(Dune::index_constant<i> domainId) const
404 { return gridGeometry(domainId).gridView(); }
405
407 template<std::size_t i>
408 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
409 { return *std::get<domainId>(gridVariablesTuple_); }
410
412 template<std::size_t i>
413 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
414 { return *std::get<domainId>(gridVariablesTuple_); }
415
418 { return *couplingManager_; }
419
422 { return *jacobian_; }
423
426 { return *residual_; }
427
429 const SolutionVector& prevSol() const
430 { return *prevSol_; }
431
436 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
437 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
438
444 { prevSol_ = &u; }
445
450 { return isStationaryProblem_; }
451
455 template<std::size_t i>
456 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
457 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
458
459protected:
461 std::shared_ptr<CouplingManager> couplingManager_;
462
463private:
464 // reset the residual vector to 0.0
465 void resetResidual_()
466 {
467 if(!residual_)
468 {
469 residual_ = std::make_shared<SolutionVector>();
470 setResidualSize(*residual_);
471 }
472
473 (*residual_) = 0.0;
474 }
475
476 // reset the jacobian vector to 0.0
477 void resetJacobian_()
478 {
479 if(!jacobian_)
480 {
481 jacobian_ = std::make_shared<JacobianMatrix>();
482 setJacobianBuildMode(*jacobian_);
483 setJacobianPattern(*jacobian_);
484 }
485
486 (*jacobian_) = 0.0;
487 }
488
489 // check if the assembler is in a correct state for assembly
490 void checkAssemblerState_() const
491 {
492 if (!isStationaryProblem_ && !prevSol_)
493 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
494
495 if (isStationaryProblem_ && prevSol_)
496 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
497 << " Did you forget to set the timeLoop to make this problem instationary?");
498 }
499
500 template<std::size_t i, class JacRow, class SubRes>
501 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
502 const SolutionVector& curSol)
503 {
504 assemble_(domainId, [&](const auto& element)
505 {
506 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
507 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
508 });
509 }
510
511 template<std::size_t i, class SubRes>
512 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
513 const SolutionVector& curSol)
514 {
515 assemble_(domainId, [&](const auto& element)
516 {
517 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
518 subDomainAssembler.assembleResidual(subRes);
519 });
520 }
521
528 template<std::size_t i, class AssembleElementFunc>
529 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
530 {
531 // let the local assembler add the element contributions
532 for (const auto& element : elements(gridView(domainId)))
533 assembleElement(element);
534 }
535
536 // get diagonal block pattern
537 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
538 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
539 Dune::index_constant<j> domainJ) const
540 {
541 const auto& gg = gridGeometry(domainI);
542 auto pattern = getJacobianPattern<isImplicit()>(gg);
543 couplingManager_->extendJacobianPattern(domainI, pattern);
544 return pattern;
545 }
546
547 // get coupling block pattern
548 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
549 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
550 Dune::index_constant<j> domainJ) const
551 {
552 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
553 domainI, gridGeometry(domainI),
554 domainJ, gridGeometry(domainJ));
555 }
556
557 // build periodic contraints into the system matrix
558 template<std::size_t i, class JacRow, class Sol, class GG>
559 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Sol& res, const GG& gridGeometry, const Sol& curSol)
560 {
561 if constexpr (GG::discMethod == DiscretizationMethods::box || GG::discMethod == DiscretizationMethods::fcstaggered)
562 {
563 for (const auto& m : gridGeometry.periodicVertexMap())
564 {
565 if (m.first < m.second)
566 {
567 auto& jac = jacRow[domainI];
568
569 // add the second row to the first
570 res[m.first] += res[m.second];
571
572 // enforce the solution of the first periodic DOF to the second one
573 res[m.second] = curSol[m.second] - curSol[m.first];
574
575 const auto end = jac[m.second].end();
576 for (auto it = jac[m.second].begin(); it != end; ++it)
577 jac[m.first][it.index()] += (*it);
578
579 // enforce constraint in second row
580 for (auto it = jac[m.second].begin(); it != end; ++it)
581 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
582
583 using namespace Dune::Hybrid;
584 forEach(makeIncompleteIntegerSequence<JacRow::size(), domainI>(), [&](const auto couplingDomainId)
585 {
586 auto& jacCoupling = jacRow[couplingDomainId];
587
588 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
589 jacCoupling[m.first][it.index()] += (*it);
590
591 for (auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
592 (*it) = 0.0;
593 });
594 }
595 }
596 }
597 }
598
600 ProblemTuple problemTuple_;
601
603 GridGeometryTuple gridGeometryTuple_;
604
606 GridVariablesTuple gridVariablesTuple_;
607
609 std::shared_ptr<const TimeLoop> timeLoop_;
610
612 const SolutionVector* prevSol_ = nullptr;
613
615 bool isStationaryProblem_;
616
618 std::shared_ptr<JacobianMatrix> jacobian_;
619 std::shared_ptr<SolutionVector> residual_;
620
622 bool warningIssued_;
623};
624
625} // end namespace Dumux
626
627#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper function to generate Jacobian pattern for different discretization methods.
Utilities for template meta programming.
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
Helper function to generate Jacobian pattern for multi domain models.
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...
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)
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:42
Definition: adapt.hh:29
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:71
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
constexpr Box box
Definition: method.hh:139
constexpr FCStaggered fcstaggered
Definition: method.hh:142
Definition: common/properties.hh:74
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
Definition: parallelhelpers.hh:430
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:64
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:388
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:369
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:87
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:88
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:84
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:421
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:95
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:220
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:86
CMType CouplingManager
Definition: multidomain/fvassembler.hh:90
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:213
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:379
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:238
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:413
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:443
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:78
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:296
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:398
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:359
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:311
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:403
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:81
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:417
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:436
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:429
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:408
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:461
MDTraits Traits
Definition: multidomain/fvassembler.hh:69
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:172
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:324
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:75
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:456
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:393
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:151
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:194
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:449
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:71
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:343
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:425
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Declares all properties used in Dumux.
Manages the handling of time dependent problems.