3.2-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
41
46
48
49namespace Dumux {
50
60template<class MDTraits, class CMType, DiffMethod diffMethod, bool useImplicitAssembly = true>
62{
63 template<std::size_t id>
64 using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
65
66public:
67 using Traits = MDTraits;
68
69 using Scalar = typename MDTraits::Scalar;
70
72 template<std::size_t id>
74
75 template<std::size_t id>
76 using GridVariables = typename MDTraits::template SubDomain<id>::GridVariables;
77
78 template<std::size_t id>
79 using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
80
81 template<std::size_t id>
82 using Problem = typename MDTraits::template SubDomain<id>::Problem;
83
84 using JacobianMatrix = typename MDTraits::JacobianMatrix;
85 using SolutionVector = typename MDTraits::SolutionVector;
87
88 using CouplingManager = CMType;
89
93 static constexpr bool isImplicit()
94 { return useImplicitAssembly; }
95
96private:
97
98 using ProblemTuple = typename MDTraits::template TupleOfSharedPtrConst<Problem>;
99 using GridGeometryTuple = typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
100 using GridVariablesTuple = typename MDTraits::template TupleOfSharedPtr<GridVariables>;
101
103 using ThisType = MultiDomainFVAssembler<MDTraits, CouplingManager, diffMethod, isImplicit()>;
104
105 template<DiscretizationMethod discMethod, std::size_t id>
106 struct SubDomainAssemblerType;
107
108 template<std::size_t id>
109 struct SubDomainAssemblerType<DiscretizationMethod::cctpfa, id>
110 {
111 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
112 };
113
114 template<std::size_t id>
115 struct SubDomainAssemblerType<DiscretizationMethod::ccmpfa, id>
116 {
117 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
118 };
119
120 template<std::size_t id>
121 struct SubDomainAssemblerType<DiscretizationMethod::box, id>
122 {
123 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
124 };
125
126 template<std::size_t id>
127 struct SubDomainAssemblerType<DiscretizationMethod::staggered, id>
128 {
129 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod, isImplicit()>;
130 };
131
132 template<std::size_t id>
133 using SubDomainAssembler = typename SubDomainAssemblerType<GridGeometry<id>::discMethod, id>::type;
134
135public:
136
137
144 GridGeometryTuple&& gridGeometry,
145 GridVariablesTuple&& gridVariables,
146 std::shared_ptr<CouplingManager> couplingManager)
148 , problemTuple_(problem)
149 , gridGeometryTuple_(gridGeometry)
150 , gridVariablesTuple_(gridVariables)
151 , timeLoop_()
152 , isStationaryProblem_(true)
153 {
154 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
155 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
156 }
157
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
200 template<class DeprecatedGridGeometryTuple, class DeprecatedGridVariablesTuple,
201 bool isStaggered = std::tuple_element_t<0, GridGeometryTuple>::element_type::discMethod == DiscretizationMethod::staggered,
202 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
203 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
204 [[deprecated("Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
206 DeprecatedGridGeometryTuple&& gridGeometry,
207 DeprecatedGridVariablesTuple&& gridVariables,
208 std::shared_ptr<CouplingManager> couplingManager)
210 makeIndexSequenceWithOffset<2,std::tuple_size<DeprecatedGridGeometryTuple>::value - 2>())
211 {}
212
213 template<class DeprecatedGridGeometryTuple, class DeprecatedGridVariablesTuple,
214 bool isStaggered = std::tuple_element_t<0, GridGeometryTuple>::element_type::discMethod == DiscretizationMethod::staggered,
215 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
216 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
217 [[deprecated("Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
219 DeprecatedGridGeometryTuple&& gridGeometry,
220 DeprecatedGridVariablesTuple&& gridVariables,
221 std::shared_ptr<CouplingManager> couplingManager,
222 std::shared_ptr<const TimeLoop> timeLoop,
223 const SolutionVector& prevSol)
224 : MultiDomainFVAssembler(std::forward<ProblemTuple>(problem), gridGeometry, gridVariables, couplingManager, timeLoop, prevSol,
225 makeIndexSequenceWithOffset<2,std::tuple_size<DeprecatedGridGeometryTuple>::value - 2>())
226 {}
227
228 template<class DeprecatedGridGeometryTuple, class DeprecatedGridVariablesTuple,
229 bool isStaggered = std::tuple_element_t<0, GridGeometryTuple>::element_type::discMethod == DiscretizationMethod::staggered,
230 bool isDeprecated = !std::is_convertible<DeprecatedGridVariablesTuple, GridVariablesTuple>::value,
231 std::enable_if_t<isStaggered && isDeprecated, int> = 0>
232 [[deprecated("Please change the order within your tuples for gridGeometry and gridVariables. Will be removed after release 3.2!")]]
234 DeprecatedGridGeometryTuple&& gridGeometry,
235 DeprecatedGridVariablesTuple&& gridVariables,
236 std::shared_ptr<CouplingManager> couplingManager,
237 std::shared_ptr<const TimeLoop> timeLoop)
238 : MultiDomainFVAssembler(std::forward<ProblemTuple>(problem), gridGeometry, gridVariables, couplingManager, timeLoop, nullptr,
239 makeIndexSequenceWithOffset<2,std::tuple_size<DeprecatedGridGeometryTuple>::value - 2>())
240 {}
241
247 {
248 checkAssemblerState_();
249 resetJacobian_();
250 resetResidual_();
251
252 using namespace Dune::Hybrid;
253 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainId)
254 {
255 auto& jacRow = (*jacobian_)[domainId];
256 auto& subRes = (*residual_)[domainId];
257 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
258 });
259 }
260
263 {
264 resetResidual_();
265 assembleResidual(*residual_, curSol);
266 }
267
270 {
271 checkAssemblerState_();
272
273 // update the grid variables for the case of active caching
274 updateGridVariables(curSol);
275
276 using namespace Dune::Hybrid;
277 forEach(integralRange(Dune::Hybrid::size(r)), [&](const auto domainId)
278 {
279 auto& subRes = r[domainId];
280 this->assembleResidual_(domainId, subRes, curSol);
281 });
282 }
283
287 {
290 assembleResidual(residual, curSol);
291
292 // calculate the square norm of the residual
293 const Scalar result2 = residual.two_norm2();
294 using std::sqrt;
295 return sqrt(result2);
296 }
297
303 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
304 std::shared_ptr<SolutionVector> r)
305 {
306 jacobian_ = A;
307 residual_ = r;
308
309 setJacobianBuildMode(*jacobian_);
310 setJacobianPattern(*jacobian_);
311 setResidualSize(*residual_);
312 }
313
319 {
320 jacobian_ = std::make_shared<JacobianMatrix>();
321 residual_ = std::make_shared<SolutionVector>();
322
323 setJacobianBuildMode(*jacobian_);
324 setJacobianPattern(*jacobian_);
325 setResidualSize(*residual_);
326 }
327
332 {
333 using namespace Dune::Hybrid;
334 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto i)
335 {
336 forEach(jac[i], [&](auto& jacBlock)
337 {
338 using BlockType = std::decay_t<decltype(jacBlock)>;
339 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
340 jacBlock.setBuildMode(BlockType::BuildMode::random);
341 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
342 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
343 });
344 });
345 }
346
351 {
352 using namespace Dune::Hybrid;
353 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](const auto domainI)
354 {
355 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](const auto domainJ)
356 {
357 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
358 pattern.exportIdx(jac[domainI][domainJ]);
359 });
360 });
361 }
362
367 {
368 using namespace Dune::Hybrid;
369 forEach(integralRange(Dune::Hybrid::size(res)), [&](const auto domainId)
370 { res[domainId].resize(this->numDofs(domainId)); });
371 }
372
377 {
378 using namespace Dune::Hybrid;
379 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
380 { this->gridVariables(domainId).update(curSol[domainId]); });
381 }
382
386 void resetTimeStep(const SolutionVector& curSol)
387 {
388 using namespace Dune::Hybrid;
389 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](const auto domainId)
390 { this->gridVariables(domainId).resetTimeStep(curSol[domainId]); });
391 }
392
394 template<std::size_t i>
395 std::size_t numDofs(Dune::index_constant<i> domainId) const
396 { return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
397
399 template<std::size_t i>
400 const auto& problem(Dune::index_constant<i> domainId) const
401 { return *std::get<domainId>(problemTuple_); }
402
404 template<std::size_t i>
405 const auto& gridGeometry(Dune::index_constant<i> domainId) const
406 { return *std::get<domainId>(gridGeometryTuple_); }
407
409 template<std::size_t i>
410 const auto& gridView(Dune::index_constant<i> domainId) const
411 { return gridGeometry(domainId).gridView(); }
412
414 template<std::size_t i>
415 GridVariables<i>& gridVariables(Dune::index_constant<i> domainId)
416 { return *std::get<domainId>(gridVariablesTuple_); }
417
419 template<std::size_t i>
420 const GridVariables<i>& gridVariables(Dune::index_constant<i> domainId) const
421 { return *std::get<domainId>(gridVariablesTuple_); }
422
425 { return *couplingManager_; }
426
429 { return *jacobian_; }
430
433 { return *residual_; }
434
436 const SolutionVector& prevSol() const
437 { return *prevSol_; }
438
443 void setTimeManager(std::shared_ptr<const TimeLoop> timeLoop)
444 { timeLoop_ = timeLoop; isStationaryProblem_ = !(static_cast<bool>(timeLoop)); }
445
451 { prevSol_ = &u; }
452
457 { return isStationaryProblem_; }
458
462 template<std::size_t i>
463 LocalResidual<i> localResidual(Dune::index_constant<i> domainId) const
464 { return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
465
466protected:
468 std::shared_ptr<CouplingManager> couplingManager_;
469
470private:
471 // reset the residual vector to 0.0
472 void resetResidual_()
473 {
474 if(!residual_)
475 {
476 residual_ = std::make_shared<SolutionVector>();
477 setResidualSize(*residual_);
478 }
479
480 (*residual_) = 0.0;
481 }
482
483 // reset the jacobian vector to 0.0
484 void resetJacobian_()
485 {
486 if(!jacobian_)
487 {
488 jacobian_ = std::make_shared<JacobianMatrix>();
489 setJacobianBuildMode(*jacobian_);
490 setJacobianPattern(*jacobian_);
491 }
492
493 (*jacobian_) = 0.0;
494 }
495
496 // check if the assembler is in a correct state for assembly
497 void checkAssemblerState_() const
498 {
499 if (!isStationaryProblem_ && !prevSol_)
500 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
501
502 if (isStationaryProblem_ && prevSol_)
503 DUNE_THROW(Dune::InvalidStateException, "Assembling stationary problem but a previous solution was set."
504 << " Did you forget to set the timeLoop to make this problem instationary?");
505 }
506
507 template<std::size_t i, class JacRow, class SubRes>
508 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
509 const SolutionVector& curSol)
510 {
511 assemble_(domainId, [&](const auto& element)
512 {
513 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
514 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
515 });
516 }
517
518 template<std::size_t i, class SubRes>
519 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
520 const SolutionVector& curSol)
521 {
522 assemble_(domainId, [&](const auto& element)
523 {
524 SubDomainAssembler<i> subDomainAssembler(*this, element, curSol, *couplingManager_);
525 subDomainAssembler.assembleResidual(subRes);
526 });
527 }
528
535 template<std::size_t i, class AssembleElementFunc>
536 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement) const
537 {
538 // let the local assembler add the element contributions
539 for (const auto& element : elements(gridView(domainId)))
540 assembleElement(element);
541 }
542
543 // get diagonal block pattern
544 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i==j), int> = 0>
545 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
546 Dune::index_constant<j> domainJ) const
547 {
548 const auto& gg = gridGeometry(domainI);
549 auto pattern = getJacobianPattern<isImplicit()>(gg);
550 couplingManager_->extendJacobianPattern(domainI, pattern);
551 return pattern;
552 }
553
554 // get coupling block pattern
555 template<std::size_t i, std::size_t j, typename std::enable_if_t<(i!=j), int> = 0>
556 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
557 Dune::index_constant<j> domainJ) const
558 {
559 return getCouplingJacobianPattern<isImplicit()>(*couplingManager_,
560 domainI, gridGeometry(domainI),
561 domainJ, gridGeometry(domainJ));
562 }
563
564 template<class DeprecatedGridGeometryTuple, class DeprecatedGridVariablesTuple, std::size_t... I>
565 [[deprecated("\n\n Please change the order within your tuples for gridGeometry and gridVariables. Usage:\n\n"
566 "auto assembler = std::make_shared<Assembler>(std::make_tuple(ffProblem, ffProblem, otherProblem, ...),\n"
567 " std::make_tuple(ffGridGeometry->faceFVGridGeometryPtr(),\n"
568 " ffFvGridGeometry->cellCenterFVGridGeometryPtr(),\n"
569 " otherFvGridGeometry, ...),\n"
570 " std::make_tuple(ffGridVariables->faceGridVariablesPtr(),\n"
571 " ffGridVariables->cellCenterGridVariablesPtr(),\n"
572 " otherGridVariables, ...),\n"
573 " couplingManager);\n\n"
574 "Will be removed after release 3.2!")]]
575 MultiDomainFVAssembler(ProblemTuple&& problem,
576 DeprecatedGridGeometryTuple&& gridGeometry,
577 DeprecatedGridVariablesTuple&& gridVariables,
578 std::shared_ptr<CouplingManager> couplingManager,
579 std::index_sequence<I...>)
581 , problemTuple_(problem)
582 , gridGeometryTuple_(std::make_tuple(std::move(std::get<1>(gridGeometry)),
583 std::move(std::get<0>(gridGeometry)),
584 std::move(std::get<I>(gridGeometry))...))
585 , gridVariablesTuple_(std::make_tuple(std::move(std::get<1>(gridVariables)),
586 std::move(std::get<0>(gridVariables)),
587 std::move(std::get<I>(gridVariables))...))
588 , timeLoop_()
589 , isStationaryProblem_(true)
590 {
591 static_assert(isImplicit(), "Explicit assembler for stationary problem doesn't make sense!");
592 std::cout << "Instantiated assembler for a stationary problem." << std::endl;
593 }
594
595 template<class DeprecatedGridGeometryTuple, class DeprecatedGridVariablesTuple, std::size_t... I>
596 [[deprecated("\n\n Please change the order within your tuples for gridGeometry and gridVariables. Usage:\n\n"
597 "auto assembler = std::make_shared<Assembler>(std::make_tuple(ffProblem, ffProblem, otherProblem, ...),\n"
598 " std::make_tuple(ffGridGeometry->faceFVGridGeometryPtr(),\n"
599 " ffFvGridGeometry->cellCenterFVGridGeometryPtr(),\n"
600 " otherFvGridGeometry, ...),\n"
601 " std::make_tuple(ffGridVariables->faceGridVariablesPtr(),\n"
602 " ffGridVariables->cellCenterGridVariablesPtr(),\n"
603 " otherGridVariables, ...),\n"
604 " couplingManager,\n"
605 " timeLoop, solOld);\n\n"
606 "Will be removed after release 3.2!")]]
607 MultiDomainFVAssembler(ProblemTuple&& problem,
608 DeprecatedGridGeometryTuple&& gridGeometry,
609 DeprecatedGridVariablesTuple&& gridVariables,
610 std::shared_ptr<CouplingManager> couplingManager,
611 std::shared_ptr<const TimeLoop> timeLoop,
612 const SolutionVector& prevSol,
613 std::index_sequence<I...>)
615 , problemTuple_(problem)
616 , gridGeometryTuple_(std::make_tuple(std::move(std::get<1>(gridGeometry)),
617 std::move(std::get<0>(gridGeometry)),
618 std::move(std::get<I>(gridGeometry))...))
619 , gridVariablesTuple_(std::make_tuple(std::move(std::get<1>(gridVariables)),
620 std::move(std::get<0>(gridVariables)),
621 std::move(std::get<I>(gridVariables))...))
622 , timeLoop_(timeLoop)
623 , prevSol_(&prevSol)
624 , isStationaryProblem_(false)
625 {
626 std::cout << "Instantiated assembler for an instationary problem." << std::endl;
627 }
628
630 ProblemTuple problemTuple_;
631
633 GridGeometryTuple gridGeometryTuple_;
634
636 GridVariablesTuple gridVariablesTuple_;
637
639 std::shared_ptr<const TimeLoop> timeLoop_;
640
642 const SolutionVector* prevSol_ = nullptr;
643
645 bool isStationaryProblem_;
646
648 std::shared_ptr<JacobianMatrix> jacobian_;
649 std::shared_ptr<SolutionVector> residual_;
650};
651
652} // end namespace Dumux
653
654#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.
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
Definition: adapt.hh:29
constexpr auto makeIndexSequenceWithOffset()
Definition: utility.hh:96
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: common/properties.hh:89
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
The default time loop for instationary simulations.
Definition: timeloop.hh:117
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:62
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:395
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
Definition: multidomain/fvassembler.hh:218
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:376
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:85
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:86
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:82
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:428
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:269
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:84
CMType CouplingManager
Definition: multidomain/fvassembler.hh:88
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:262
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:386
Scalar residualNorm(const SolutionVector &curSol)
Definition: multidomain/fvassembler.hh:286
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop)
Definition: multidomain/fvassembler.hh:233
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:420
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:450
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:76
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< SolutionVector > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: multidomain/fvassembler.hh:303
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:405
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:366
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:318
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:410
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 MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:79
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
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:424
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:443
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:436
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:415
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:468
MDTraits Traits
Definition: multidomain/fvassembler.hh:67
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:331
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:73
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:143
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:463
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:400
MultiDomainFVAssembler(ProblemTuple &&problem, DeprecatedGridGeometryTuple &&gridGeometry, DeprecatedGridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
Definition: multidomain/fvassembler.hh:205
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:246
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:456
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:69
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:350
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:432
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:268
Declares all properties used in Dumux.