3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomainstaggeredlocalassembler.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_STAGGERED_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_STAGGERED_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/reservedvector.hh>
30#include <dune/common/indices.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh> // for GhostEntity
33
42
43namespace Dumux {
44
56template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool isImplicit = true>
57class SubDomainStaggeredLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler,Implementation, isImplicit>
58{
60
63 using SolutionVector = typename Assembler::SolutionVector;
64
66 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
67 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
68 using Scalar = typename GridVariables::Scalar;
69
70 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
71 using CellCenterResidualValue = typename LocalResidual::CellCenterResidualValue;
72 using FaceResidualValue = typename LocalResidual::FaceResidualValue;
73
74 using GridGeometry = typename GridVariables::GridGeometry;
75 using FVElementGeometry = typename GridGeometry::LocalView;
76 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
77 using GridView = typename GridGeometry::GridView;
78 using Element = typename GridView::template Codim<0>::Entity;
79
80 using CouplingManager = typename Assembler::CouplingManager;
81
82 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
83
84public:
85 static constexpr auto domainId = typename Dune::index_constant<id>();
86 static constexpr auto cellCenterId = GridGeometry::cellCenterIdx();
87 static constexpr auto faceId = GridGeometry::faceIdx();
88
89 static constexpr auto numEqCellCenter = CellCenterResidualValue::dimension;
90 static constexpr auto faceOffset = numEqCellCenter;
91
92 using ParentType::ParentType;
93
95 const Element& element,
96 const SolutionVector& curSol,
97 CouplingManager& couplingManager)
99 element,
100 curSol,
101 localView(assembler.gridGeometry(domainId)),
102 localView(assembler.gridVariables(domainId).curGridVolVars()),
103 localView(assembler.gridVariables(domainId).prevGridVolVars()),
104 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
106 (element.partitionType() == Dune::GhostEntity))
107 , curElemFaceVars_(localView(assembler.gridVariables(domainId).curGridFaceVars()))
108 , prevElemFaceVars_(localView(assembler.gridVariables(domainId).prevGridFaceVars()))
109 , couplingManager_(couplingManager)
110 {}
111
116 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
117 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
118 {
119 this->asImp_().bindLocalViews();
120 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
121
122 assembleJacobianAndResidualImpl_(domainId, jacRow, res, gridVariables);
123 }
124
128 template<class SubSol>
129 void assembleResidual(SubSol& res)
130 {
131 this->asImp_().bindLocalViews();
132 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
133
134 if constexpr (domainId == cellCenterId)
135 {
136 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
137 res[cellCenterGlobalI] = this->asImp_().assembleCellCenterResidualImpl();
138 }
139 else
140 {
141 for (auto&& scvf : scvfs(this->fvGeometry()))
142 res[scvf.dofIndex()] += this->asImp_().assembleFaceResidualImpl(scvf);
143 }
144 }
145
150 CellCenterResidualValue evalLocalResidualForCellCenter() const
151 {
152 if (!isImplicit)
153 if (this->assembler().isStationaryProblem())
154 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
155
156 if (this->elementIsGhost())
157 {
158 return CellCenterResidualValue(0.0);
159 }
160
163 }
164
170 CellCenterResidualValue evalLocalResidualForCellCenter(const ElementVolumeVariables& elemVolVars,
171 const ElementFaceVariables& elemFaceVars) const
172 {
173 auto residual = evalLocalFluxAndSourceResidualForCellCenter(elemVolVars, elemFaceVars);
174
175 if (!this->assembler().isStationaryProblem())
177
178 // handle cells with a fixed Dirichlet value
179 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
180 const auto& scvI = this->fvGeometry().scv(cellCenterGlobalI);
181 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
182 {
183 static constexpr auto offset = numEq - numEqCellCenter;
184 if (this->problem().isDirichletCell(this->element(), this->fvGeometry(), scvI, pvIdx + offset))
185 residual[pvIdx] = this->curSol()[cellCenterId][cellCenterGlobalI][pvIdx] - this->problem().dirichlet(this->element(), scvI)[pvIdx + offset];
186 }
187
188 return residual;
189 }
190
196 CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter() const
197 {
200 }
201
209 CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter(const ElementVolumeVariables& elemVolVars, const ElementFaceVariables& elemFaceVars) const
210 {
211 return this->localResidual().evalFluxAndSourceForCellCenter(this->element(), this->fvGeometry(), elemVolVars, elemFaceVars, this->elemBcTypes(), this->elemFluxVarsCache());
212 }
213
219 CellCenterResidualValue evalLocalStorageResidualForCellCenter() const
220 {
221 return this->localResidual().evalStorageForCellCenter(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars());
222 }
223
229 FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace& scvf) const
230 {
231 if (!isImplicit)
232 if (this->assembler().isStationaryProblem())
233 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
234
235 if (this->elementIsGhost())
236 {
237 return FaceResidualValue(0.0);
238 }
239
242 }
243
250 FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace& scvf,
251 const ElementVolumeVariables& elemVolVars,
252 const ElementFaceVariables& elemFaceVars) const
253 {
254 auto residual = evalLocalFluxAndSourceResidualForFace(scvf, elemVolVars, elemFaceVars);
255
256 if (!this->assembler().isStationaryProblem())
257 residual += evalLocalStorageResidualForFace(scvf);
258
259 this->localResidual().evalDirichletBoundariesForFace(residual, this->problem(), this->element(),
260 this->fvGeometry(), scvf, elemVolVars, elemFaceVars,
261 this->elemBcTypes(), this->elemFluxVarsCache());
262
263 return residual;
264 }
265
272 FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace& scvf) const
273 {
276 }
277
285 FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace& scvf,
286 const ElementVolumeVariables& elemVolVars,
287 const ElementFaceVariables& elemFaceVars) const
288 {
289 return this->localResidual().evalFluxAndSourceForFace(this->element(), this->fvGeometry(), elemVolVars, elemFaceVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf);
290 }
291
298 FaceResidualValue evalLocalStorageResidualForFace(const SubControlVolumeFace& scvf) const
299 {
300 return this->localResidual().evalStorageForFace(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars(), this->prevElemFaceVars(), this->curElemFaceVars(), scvf);
301 }
302
303 const Problem& problem() const
304 { return this->assembler().problem(domainId); }
305
307 ElementFaceVariables& curElemFaceVars()
308 { return curElemFaceVars_; }
309
311 ElementFaceVariables& prevElemFaceVars()
312 { return prevElemFaceVars_; }
313
315 const ElementFaceVariables& curElemFaceVars() const
316 { return curElemFaceVars_; }
317
319 const ElementFaceVariables& prevElemFaceVars() const
320 { return prevElemFaceVars_; }
321
322 CouplingManager& couplingManager()
323 { return couplingManager_; }
324
325private:
326
328 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
329 auto assembleJacobianAndResidualImpl_(Dune::index_constant<cellCenterId>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
330 {
331 auto& gridVariablesI = *std::get<domainId>(gridVariables);
332 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
333 const auto residual = this->asImp_().assembleCellCenterJacobianAndResidualImpl(jacRow[domainId], gridVariablesI);
334 res[cellCenterGlobalI] = residual;
335
336
337 // for the coupling blocks
338 using namespace Dune::Hybrid;
339 static constexpr auto otherDomainIds = makeIncompleteIntegerSequence<JacobianMatrixRow::size(), domainId>{};
340 forEach(otherDomainIds, [&](auto&& domainJ)
341 {
342 this->asImp_().assembleJacobianCellCenterCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
343 });
344
345 // handle cells with a fixed Dirichlet value
346 incorporateDirichletCells_(jacRow);
347 }
348
350 template<class JacobianMatrixRow, class SubSol, class GridVariablesTuple>
351 void assembleJacobianAndResidualImpl_(Dune::index_constant<faceId>, JacobianMatrixRow& jacRow, SubSol& res, GridVariablesTuple& gridVariables)
352 {
353 auto& gridVariablesI = *std::get<domainId>(gridVariables);
354 const auto residual = this->asImp_().assembleFaceJacobianAndResidualImpl(jacRow[domainId], gridVariablesI);
355
356 for(auto&& scvf : scvfs(this->fvGeometry()))
357 res[scvf.dofIndex()] += residual[scvf.localFaceIdx()];
358
359 // for the coupling blocks
360 using namespace Dune::Hybrid;
361 static constexpr auto otherDomainIds = makeIncompleteIntegerSequence<JacobianMatrixRow::size(), domainId>{};
362 forEach(otherDomainIds, [&](auto&& domainJ)
363 {
364 this->asImp_().assembleJacobianFaceCoupling(domainJ, jacRow[domainJ], residual, gridVariablesI);
365 });
366 }
367
369 template<class JacobianMatrixRow>
370 void incorporateDirichletCells_(JacobianMatrixRow& jacRow)
371 {
372 const auto cellCenterGlobalI = problem().gridGeometry().elementMapper().index(this->element());
373
374 // overwrite the partial derivative with zero in case a fixed Dirichlet BC is used
375 static constexpr auto offset = numEq - numEqCellCenter;
376 for (int eqIdx = 0; eqIdx < numEqCellCenter; ++eqIdx)
377 {
378 if (problem().isDirichletCell(this->element(), this->fvGeometry(), this->fvGeometry().scv(cellCenterGlobalI), eqIdx + offset))
379 {
380 using namespace Dune::Hybrid;
381 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&, domainId = domainId](auto&& i)
382 {
383 auto& ccRowI = jacRow[i][cellCenterGlobalI];
384 for (auto col = ccRowI.begin(); col != ccRowI.end(); ++col)
385 {
386 ccRowI[col.index()][eqIdx] = 0.0;
387 // set the diagonal entry to 1.0
388 if ((i == domainId) && (col.index() == cellCenterGlobalI))
389 ccRowI[col.index()][eqIdx][eqIdx] = 1.0;
390 }
391 });
392 }
393 }
394 }
395
396 ElementFaceVariables curElemFaceVars_;
397 ElementFaceVariables prevElemFaceVars_;
398 CouplingManager& couplingManager_;
399};
400
411template<std::size_t id, class TypeTag, class Assembler, class Implementation>
412class SubDomainStaggeredLocalAssemblerImplicitBase : public SubDomainStaggeredLocalAssemblerBase<id, TypeTag, Assembler, Implementation>
413{
415 static constexpr auto domainId = Dune::index_constant<id>();
416public:
417 using ParentType::ParentType;
418
420 {
421 // get some references for convenience
422 auto& couplingManager = this->couplingManager();
423 const auto& element = this->element();
424 const auto& curSol = this->curSol();
425 auto&& fvGeometry = this->fvGeometry();
426 auto&& curElemVolVars = this->curElemVolVars();
427 auto&& curElemFaceVars = this->curElemFaceVars();
428 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
429
430 // bind the caches
431 couplingManager.bindCouplingContext(domainId, element, this->assembler());
432 fvGeometry.bind(element);
436 if (!this->assembler().isStationaryProblem())
437 {
438 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol());
439 this->prevElemFaceVars().bindElement(element, fvGeometry, this->assembler().prevSol());
440 }
441 }
442};
443
455template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
457
478template<std::size_t id, class TypeTag, class Assembler>
479class SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
480: public SubDomainStaggeredLocalAssemblerImplicitBase<id, TypeTag, Assembler,
481 SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, true> >
482{
483 using ThisType = SubDomainStaggeredLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
487 using CellCenterResidualValue = typename LocalResidual::CellCenterResidualValue;
488 using FaceResidualValue = typename LocalResidual::FaceResidualValue;
489 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
491 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
492 using FaceVariables = typename ElementFaceVariables::FaceVariables;
494 using FVElementGeometry = typename GridGeometry::LocalView;
495 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
497 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
500
502 static constexpr int maxNeighbors = 4*(2*ModelTraits::dim());
503 static constexpr auto domainI = Dune::index_constant<id>();
504 static constexpr auto cellCenterId = GridGeometry::cellCenterIdx();
505 static constexpr auto faceId = GridGeometry::faceIdx();
506
507 static constexpr auto numEq = ModelTraits::numEq();
508 static constexpr auto numEqCellCenter = CellCenterPrimaryVariables::dimension;
509 static constexpr auto numEqFace = FacePrimaryVariables::dimension;
510
511public:
513
514 CellCenterResidualValue assembleCellCenterResidualImpl()
515 {
516 return this->evalLocalResidualForCellCenter();
517 }
518
519 FaceResidualValue assembleFaceResidualImpl(const SubControlVolumeFace& scvf)
520 {
521 return this->evalLocalResidualForFace(scvf);
522 }
523
530 template<class JacobianMatrixDiagBlock, class GridVariables>
531 CellCenterResidualValue assembleCellCenterJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
532 {
533 assert(domainI == cellCenterId);
534
535 // get some aliases for convenience
536 const auto& element = this->element();
537 const auto& fvGeometry = this->fvGeometry();
538 auto&& curElemVolVars = this->curElemVolVars();
539 const auto& gridGeometry = this->problem().gridGeometry();
540 const auto& curSol = this->curSol()[domainI];
541
542 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
543 const auto origResidual = this->evalLocalResidualForCellCenter();
544
546 // Calculate derivatives of all cell center residuals in the element w.r.t. to other cell center dofs. //
548
549 // lambda to evaluate the derivatives for cell center dofs with respect to neighbor cells
550 auto evaluateCellCenterDerivatives = [&](const std::size_t globalJ)
551 {
552 // get the volVars of the element with respect to which we are going to build the derivative
553 auto&& scvJ = fvGeometry.scv(globalJ);
554 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
555 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ);
556 const auto origVolVars(curVolVars);
557
558 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
559 {
560 CellCenterPrimaryVariables cellCenterPriVars = curSol[globalJ];
561 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
562 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
563
564 constexpr auto offset = numEq - numEqCellCenter;
565
566 auto evalResidual = [&](Scalar priVar)
567 {
568 // update the volume variables
569 priVars[pvIdx + offset] = priVar;
570 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
571 curVolVars.update(elemSol, this->problem(), elementJ, scvJ);
572
573 // update the coupling context
574 cellCenterPriVars[pvIdx] = priVar;
575 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, cellCenterPriVars, pvIdx);
576
577 // compute element residual
578 return this->evalLocalResidualForCellCenter();
579 };
580
581 // create the vector storing the partial derivatives
582 CellCenterResidualValue partialDeriv(0.0);
583
584 // derive the residuals numerically
585 const auto& paramGroup = this->problem().paramGroup();
586 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
587 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
588 NumericDifferentiation::partialDerivative(evalResidual, priVars[pvIdx + offset], partialDeriv, origResidual,
589 eps(priVars[pvIdx + offset], pvIdx), numDiffMethod);
590
591 // update the global jacobian matrix with the current partial derivatives
592 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
593
594 // restore the original volVars
595 curVolVars = origVolVars;
596
597 // restore the undeflected state of the coupling context
598 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, curSol[globalJ], pvIdx);
599 }
600 };
601
602 // get the list of cell center dofs that have an influence on the cell center resdiual of the current element
603 const auto& connectivityMap = gridGeometry.connectivityMap();
604
605 // evaluate derivatives w.r.t. own dof
606 evaluateCellCenterDerivatives(cellCenterGlobalI);
607
608 // evaluate derivatives w.r.t. all other related cell center dofs
609 for (const auto& globalJ : connectivityMap(cellCenterId, cellCenterId, cellCenterGlobalI))
610 evaluateCellCenterDerivatives(globalJ);
611
612 return origResidual;
613 }
614
621 template<class JacobianMatrixDiagBlock, class GridVariables>
622 auto assembleFaceJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
623 {
624 assert(domainI == faceId);
625
626 // get some aliases for convenience
627 const auto& problem = this->problem();
628 const auto& element = this->element();
629 const auto& fvGeometry = this->fvGeometry();
630 const auto& gridGeometry = this->problem().gridGeometry();
631 const auto& curSol = this->curSol()[domainI];
632
633 using FaceSolutionVector = GetPropType<TypeTag, Properties::FaceSolutionVector>; // TODO: use reserved vector
634 FaceSolutionVector origResiduals;
635 origResiduals.resize(fvGeometry.numScvf());
636 origResiduals = 0.0;
637
638 // treat the local residua of the face dofs:
639 for (auto&& scvf : scvfs(fvGeometry))
640 origResiduals[scvf.localFaceIdx()] = this->evalLocalResidualForFace(scvf);
641
643 // Calculate derivatives of all face residuals in the element w.r.t. to other face dofs. //
644 // Note that we do an element-wise assembly, therefore this is only the contribution of the //
645 // current element while the contribution of the element on the opposite side of the scvf will //
646 // be added separately. //
648
649 for (auto&& scvf : scvfs(fvGeometry))
650 {
651 // set the actual dof index
652 const auto faceGlobalI = scvf.dofIndex();
653
655 const auto origFaceSolution = FaceSolution(scvf, curSol, gridGeometry);
656
657 // Lambda to evaluate the derivatives for faces
658 auto evaluateFaceDerivatives = [&](const std::size_t globalJ)
659 {
660 // get the faceVars of the face with respect to which we are going to build the derivative
661 auto& faceVars = getFaceVarAccess_(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvf);
662 const auto origFaceVars = faceVars;
663
664 for (int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
665 {
666 auto faceSolution = origFaceSolution;
667
668 auto evalResidual = [&](Scalar priVar)
669 {
670 // update the volume variables
671 faceSolution[globalJ][pvIdx] = priVar;
672 faceVars.update(faceSolution, problem, element, fvGeometry, scvf);
673
674 // update the coupling context
675 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, faceSolution[globalJ], pvIdx);
676
677 // compute face residual
678 return this->evalLocalResidualForFace(scvf);
679 };
680
681 // derive the residuals numerically
682 FaceResidualValue partialDeriv(0.0);
683 const auto& paramGroup = problem.paramGroup();
684 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
685 static const auto eps = this->couplingManager().numericEpsilon(domainI, paramGroup);
686 NumericDifferentiation::partialDerivative(evalResidual, faceSolution[globalJ][pvIdx], partialDeriv, origResiduals[scvf.localFaceIdx()],
687 eps(faceSolution[globalJ][pvIdx], pvIdx), numDiffMethod);
688
689 // update the global jacobian matrix with the current partial derivatives
690 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
691
692 // restore the original faceVars
693 faceVars = origFaceVars;
694
695 // restore the undeflected state of the coupling context
696 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalJ, origFaceSolution[globalJ], pvIdx);
697 }
698 };
699
700 // evaluate derivatives w.r.t. own dof
701 evaluateFaceDerivatives(scvf.dofIndex());
702
703 // get the list of face dofs that have an influence on the resdiual of the current face
704 const auto& connectivityMap = gridGeometry.connectivityMap();
705
706 // evaluate derivatives w.r.t. all other related face dofs
707 for (const auto& globalJ : connectivityMap(faceId, faceId, scvf.index()))
708 evaluateFaceDerivatives(globalJ);
709 }
710
711 return origResiduals;
712 }
713
718 template<class JacobianBlock, class GridVariables>
719 void assembleJacobianCellCenterCoupling(Dune::index_constant<faceId> domainJ, JacobianBlock& A,
720 const CellCenterResidualValue& origResidual, GridVariables& gridVariables)
721 {
723 // Calculate derivatives of all cell center residuals in the element w.r.t. to all coupled faces dofs //
725
726 // get some aliases for convenience
727 const auto& element = this->element();
728 const auto& fvGeometry = this->fvGeometry();
729 const auto& gridGeometry = this->problem().gridGeometry();
730 const auto& curSol = this->curSol()[domainJ];
731 // build derivatives with for cell center dofs w.r.t. cell center dofs
732 const auto cellCenterGlobalI = gridGeometry.elementMapper().index(element);
733
734 for (const auto& scvfJ : scvfs(fvGeometry))
735 {
736 const auto globalJ = scvfJ.dofIndex();
737
738 // get the faceVars of the face with respect to which we are going to build the derivative
739 auto& faceVars = getFaceVarAccess_(gridVariables.curGridFaceVars(), this->curElemFaceVars(), scvfJ);
740 const auto origFaceVars(faceVars);
741
742 for (int pvIdx = 0; pvIdx < numEqFace; ++pvIdx)
743 {
744 auto facePriVars = curSol[globalJ];
745
746 auto evalResidual = [&](Scalar priVar)
747 {
748 // update the face variables
749 facePriVars[pvIdx] = priVar;
750 faceVars.updateOwnFaceOnly(facePriVars);
751
752 // update the coupling context
753 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, facePriVars, pvIdx);
754
755 // compute element residual
756 return this->evalLocalResidualForCellCenter();
757 };
758
759 // create the vector storing the partial derivatives
760 CellCenterResidualValue partialDeriv(0.0);
761
762 // derive the residuals numerically
763 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
764 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
765 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
766 NumericDifferentiation::partialDerivative(evalResidual, facePriVars[pvIdx], partialDeriv, origResidual,
767 epsCoupl(facePriVars[pvIdx], pvIdx), numDiffMethod);
768
769 // update the global jacobian matrix with the current partial derivatives
770 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
771
772 // restore the original faceVars
773 faceVars = origFaceVars;
774
775 // restore the undeflected state of the coupling context
776 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, curSol[globalJ], pvIdx);
777 }
778 }
779 }
780
781 template<std::size_t otherId, class JacobianBlock, class GridVariables>
782 void assembleJacobianCellCenterCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
783 const CellCenterResidualValue& res, GridVariables& gridVariables)
784 {
786 // Calculate derivatives of all cell center residuals in the element w.r.t. all dofs in the coupling stencil. //
788
789 // get some aliases for convenience
790 const auto& element = this->element();
791
792 // get stencil information
793 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
794
795 if (stencil.empty())
796 return;
797
798 for (const auto globalJ : stencil)
799 {
800 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
801 const auto& curSol = this->curSol()[domainJ];
802 const auto origPriVarsJ = curSol[globalJ];
803
804 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
805 {
806 auto evalCouplingResidual = [&](Scalar priVar)
807 {
808 auto deflectedPriVarsJ = origPriVarsJ;
809 deflectedPriVarsJ[pvIdx] = priVar;
810 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedPriVarsJ, pvIdx);
811 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
812 };
813
814 // create the vector storing the partial derivatives
815 CellCenterResidualValue partialDeriv(0.0);
816
817 // derive the residuals numerically
818 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
819 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
820 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
821 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDeriv, origResidual,
822 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
823
824 // update the global stiffness matrix with the current partial derivatives
825 const auto cellCenterGlobalI = this->problem().gridGeometry().elementMapper().index(element);
826 updateGlobalJacobian_(A, cellCenterGlobalI, globalJ, pvIdx, partialDeriv);
827
828 // restore the undeflected state of the coupling context
829 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origPriVarsJ, pvIdx);
830 }
831 }
832 }
833
838 template<class JacobianBlock, class ElementResidualVector, class GridVariables>
839 void assembleJacobianFaceCoupling(Dune::index_constant<cellCenterId> domainJ, JacobianBlock& A,
840 const ElementResidualVector& origResiduals, GridVariables& gridVariables)
841 {
843 // Calculate derivatives of all face residuals in the element w.r.t. all coupled cell center dofs //
845
846 // get some aliases for convenience
847 const auto& problem = this->problem();
848 const auto& fvGeometry = this->fvGeometry();
849 const auto& gridGeometry = this->problem().gridGeometry();
850 const auto& connectivityMap = gridGeometry.connectivityMap();
851 const auto& curSol = this->curSol()[domainJ];
852
853 // build derivatives with for cell center dofs w.r.t. cell center dofs
854 for (auto&& scvf : scvfs(fvGeometry))
855 {
856 // set the actual dof index
857 const auto faceGlobalI = scvf.dofIndex();
858
859 // build derivatives with for face dofs w.r.t. cell center dofs
860 for (const auto& globalJ : connectivityMap(faceId, cellCenterId, scvf.index()))
861 {
862 // get the volVars of the element with respect to which we are going to build the derivative
863 auto&& scvJ = fvGeometry.scv(globalJ);
864 const auto elementJ = fvGeometry.gridGeometry().element(globalJ);
865 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), this->curElemVolVars(), scvJ);
866 const auto origVolVars(curVolVars);
867 const auto origCellCenterPriVars = curSol[globalJ];
868
869 for (int pvIdx = 0; pvIdx < numEqCellCenter; ++pvIdx)
870 {
871 using PrimaryVariables = typename VolumeVariables::PrimaryVariables;
872 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(origCellCenterPriVars);
873
874 constexpr auto offset = PrimaryVariables::dimension - CellCenterPrimaryVariables::dimension;
875
876 auto evalResidual = [&](Scalar priVar)
877 {
878 // update the volume variables
879 priVars[pvIdx + offset] = priVar;
880 auto elemSol = elementSolution<FVElementGeometry>(std::move(priVars));
881 curVolVars.update(elemSol, problem, elementJ, scvJ);
882
883 // update the coupling context
884 auto deflectedCellCenterPriVars = origCellCenterPriVars;
885 deflectedCellCenterPriVars[pvIdx] = priVar;
886 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedCellCenterPriVars, pvIdx);
887
888 // compute face residual
889 return this->evalLocalResidualForFace(scvf);
890 };
891
892 // derive the residuals numerically
893 FaceResidualValue partialDeriv(0.0);
894 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
895 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
896 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
897 NumericDifferentiation::partialDerivative(evalResidual, priVars[pvIdx + offset], partialDeriv, origResiduals[scvf.localFaceIdx()],
898 epsCoupl(priVars[pvIdx + offset], pvIdx), numDiffMethod);
899
900 // update the global jacobian matrix with the current partial derivatives
901 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
902
903 // restore the original volVars
904 curVolVars = origVolVars;
905
906 // restore the undeflected state of the coupling context
907 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origCellCenterPriVars, pvIdx);
908 }
909 }
910 }
911 }
912
913 template<std::size_t otherId, class JacobianBlock, class ElementResidualVector, class GridVariables>
914 void assembleJacobianFaceCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
915 const ElementResidualVector& res, GridVariables& gridVariables)
916 {
918 // Calculate derivatives of all face residuals in the element w.r.t. all dofs in the coupling stencil. //
920
921 // get some aliases for convenience
922 const auto& fvGeometry = this->fvGeometry();
923 const auto& curSol = this->curSol()[domainJ];
924
925 // build derivatives with for cell center dofs w.r.t. cell center dofs
926 for (auto&& scvf : scvfs(fvGeometry))
927 {
928 // set the actual dof index
929 const auto faceGlobalI = scvf.dofIndex();
930
931 // get stencil information
932 const auto& stencil = this->couplingManager().couplingStencil(domainI, scvf, domainJ);
933
934 if (stencil.empty())
935 continue;
936
937 // build derivatives with for face dofs w.r.t. cell center dofs
938 for (const auto& globalJ : stencil)
939 {
940 const auto origPriVarsJ = curSol[globalJ];
941 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, scvf, *this, domainJ, globalJ);
942
943 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
944 {
945 auto evalCouplingResidual = [&](Scalar priVar)
946 {
947 auto deflectedPriVars = origPriVarsJ;
948 deflectedPriVars[pvIdx] = priVar;
949 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, deflectedPriVars, pvIdx);
950 return this->couplingManager().evalCouplingResidual(domainI, scvf, *this, domainJ, globalJ);
951 };
952
953 // derive the residuals numerically
954 FaceResidualValue partialDeriv(0.0);
955 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
956 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
957 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
958 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDeriv, origResidual,
959 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
960
961 // update the global stiffness matrix with the current partial derivatives
962 updateGlobalJacobian_(A, faceGlobalI, globalJ, pvIdx, partialDeriv);
963
964 // restore the undeflected state of the coupling context
965 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, origPriVarsJ, pvIdx);
966 }
967 }
968 }
969 }
970
971 template<class JacobianMatrixDiagBlock, class GridVariables>
972 void evalAdditionalDerivatives(const std::vector<std::size_t>& additionalDofDependencies,
973 JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
974 { }
975
981 template<class SubMatrix, class CCOrFacePrimaryVariables>
982 static void updateGlobalJacobian_(SubMatrix& matrix,
983 const int globalI,
984 const int globalJ,
985 const int pvIdx,
986 const CCOrFacePrimaryVariables& partialDeriv)
987 {
988 for (int eqIdx = 0; eqIdx < partialDeriv.size(); eqIdx++)
989 {
990 // A[i][col][eqIdx][pvIdx] is the rate of change of
991 // the residual of equation 'eqIdx' at dof 'i'
992 // depending on the primary variable 'pvIdx' at dof
993 // 'col'.
994
995 assert(pvIdx >= 0);
996 assert(eqIdx < matrix[globalI][globalJ].size());
997 assert(pvIdx < matrix[globalI][globalJ][eqIdx].size());
998 matrix[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
999 }
1000 }
1001
1002private:
1003
1004 FaceVariables& getFaceVarAccess_(GridFaceVariables& gridFaceVariables, ElementFaceVariables& elemFaceVars, const SubControlVolumeFace& scvf)
1005 {
1007 return gridFaceVariables.faceVars(scvf.index());
1008 else
1009 return elemFaceVars[scvf];
1010 }
1011};
1012
1013} // end namespace Dumux
1014
1015#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
Utilities for template meta programming.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
A arithmetic block vector type based on DUNE's reserved vector.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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 GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Definition: deprecated.hh:149
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:273
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: fvlocalassemblerbase.hh:113
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:61
A base class for all multidomain local assemblers (staggered)
Definition: subdomainstaggeredlocalassembler.hh:58
FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate the local residual for the current face. Automatically chooses the t...
Definition: subdomainstaggeredlocalassembler.hh:229
const ElementFaceVariables & curElemFaceVars() const
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:315
FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: subdomainstaggeredlocalassembler.hh:285
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSol &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainstaggeredlocalassembler.hh:117
ElementFaceVariables & prevElemFaceVars()
The element volume variables of the provious time step.
Definition: subdomainstaggeredlocalassembler.hh:311
static constexpr auto domainId
Definition: subdomainstaggeredlocalassembler.hh:85
void assembleResidual(SubSol &res)
Assemble the residual only.
Definition: subdomainstaggeredlocalassembler.hh:129
const Problem & problem() const
Definition: subdomainstaggeredlocalassembler.hh:303
SubDomainStaggeredLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainstaggeredlocalassembler.hh:94
FaceResidualValue evalLocalFluxAndSourceResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: subdomainstaggeredlocalassembler.hh:272
CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter() const
Convenience function to evaluate the flux and source terms (i.e, the terms without a time derivative)...
Definition: subdomainstaggeredlocalassembler.hh:196
static constexpr auto faceId
Definition: subdomainstaggeredlocalassembler.hh:87
static constexpr auto faceOffset
Definition: subdomainstaggeredlocalassembler.hh:90
FaceResidualValue evalLocalResidualForFace(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the complete local residual for the current face.
Definition: subdomainstaggeredlocalassembler.hh:250
static constexpr auto numEqCellCenter
Definition: subdomainstaggeredlocalassembler.hh:89
const ElementFaceVariables & prevElemFaceVars() const
The element volume variables of the provious time step.
Definition: subdomainstaggeredlocalassembler.hh:319
CellCenterResidualValue evalLocalResidualForCellCenter() const
Convenience function to evaluate the complete local residual for the current element....
Definition: subdomainstaggeredlocalassembler.hh:150
FaceResidualValue evalLocalStorageResidualForFace(const SubControlVolumeFace &scvf) const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: subdomainstaggeredlocalassembler.hh:298
CouplingManager & couplingManager()
Definition: subdomainstaggeredlocalassembler.hh:322
static constexpr auto cellCenterId
Definition: subdomainstaggeredlocalassembler.hh:86
CellCenterResidualValue evalLocalResidualForCellCenter(const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the complete local residual for the current cell center.
Definition: subdomainstaggeredlocalassembler.hh:170
CellCenterResidualValue evalLocalFluxAndSourceResidualForCellCenter(const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Evaluates the flux and source terms (i.e, the terms without a time derivative) of the local residual ...
Definition: subdomainstaggeredlocalassembler.hh:209
CellCenterResidualValue evalLocalStorageResidualForCellCenter() const
Convenience function to evaluate storage term (i.e, the term with a time derivative) of the local res...
Definition: subdomainstaggeredlocalassembler.hh:219
ElementFaceVariables & curElemFaceVars()
The current element volume variables.
Definition: subdomainstaggeredlocalassembler.hh:307
A base class for all implicit multidomain local assemblers (staggered)
Definition: subdomainstaggeredlocalassembler.hh:413
void bindLocalViews()
Definition: subdomainstaggeredlocalassembler.hh:419
The staggered multidomain local assembler.
Definition: subdomainstaggeredlocalassembler.hh:456
Staggered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: subdomainstaggeredlocalassembler.hh:482
void assembleJacobianFaceCoupling(Dune::index_constant< cellCenterId > domainJ, JacobianBlock &A, const ElementResidualVector &origResiduals, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:839
void assembleJacobianFaceCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:914
FaceResidualValue assembleFaceResidualImpl(const SubControlVolumeFace &scvf)
Definition: subdomainstaggeredlocalassembler.hh:519
CellCenterResidualValue assembleCellCenterJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:531
void assembleJacobianCellCenterCoupling(Dune::index_constant< faceId > domainJ, JacobianBlock &A, const CellCenterResidualValue &origResidual, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:719
CellCenterResidualValue assembleCellCenterResidualImpl()
Definition: subdomainstaggeredlocalassembler.hh:514
static void updateGlobalJacobian_(SubMatrix &matrix, const int globalI, const int globalJ, const int pvIdx, const CCOrFacePrimaryVariables &partialDeriv)
Updates the current global Jacobian matrix with the partial derivatives of all equations in regard to...
Definition: subdomainstaggeredlocalassembler.hh:982
void evalAdditionalDerivatives(const std::vector< std::size_t > &additionalDofDependencies, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:972
auto assembleFaceJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainstaggeredlocalassembler.hh:622
void assembleJacobianCellCenterCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const CellCenterResidualValue &res, GridVariables &gridVariables)
Definition: subdomainstaggeredlocalassembler.hh:782
Declares all properties used in Dumux.
The local element solution class for staggered methods.