3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomaincclocalassembler.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_CC_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/indices.hh>
30#include <dune/common/reservedvector.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh> // for GhostEntity
33#include <dune/istl/matrixindexset.hh>
34
45
46namespace Dumux {
47
59template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool implicit>
60class SubDomainCCLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler,Implementation, implicit>
61{
63
67 using SolutionVector = typename Assembler::SolutionVector;
70
72 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
73 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
74 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
75 using Scalar = typename GridVariables::Scalar;
76
77 using GridGeometry = typename GridVariables::GridGeometry;
78 using FVElementGeometry = typename GridGeometry::LocalView;
79 using SubControlVolume = typename GridGeometry::SubControlVolume;
80 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
81 using Extrusion = Extrusion_t<GridGeometry>;
82 using GridView = typename GridGeometry::GridView;
83 using Element = typename GridView::template Codim<0>::Entity;
84
85 using CouplingManager = typename Assembler::CouplingManager;
86 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
87
88public:
90 static constexpr auto domainId = typename Dune::index_constant<id>();
94 using ParentType::ParentType;
95
97 explicit SubDomainCCLocalAssemblerBase(const Assembler& assembler,
98 const Element& element,
99 const SolutionVector& curSol,
100 CouplingManager& couplingManager)
102 element,
103 curSol,
104 localView(assembler.gridGeometry(domainId)),
105 localView(assembler.gridVariables(domainId).curGridVolVars()),
106 localView(assembler.gridVariables(domainId).prevGridVolVars()),
107 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
109 (element.partitionType() == Dune::GhostEntity))
110 , couplingManager_(couplingManager)
111 {}
112
117 template<class JacobianMatrixRow, class GridVariablesTuple>
118 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSolutionVector& res, GridVariablesTuple& gridVariables)
119 {
120 this->asImp_().bindLocalViews();
121
122 // for the diagonal jacobian block
123 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
124 res[globalI] = this->asImp_().assembleJacobianAndResidualImplInverse(jacRow[domainId], *std::get<domainId>(gridVariables));
125
126 // for the coupling blocks
127 using namespace Dune::Hybrid;
128 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
129 {
130 if (i != id)
131 this->assembleJacobianCoupling(i, jacRow, res[globalI], gridVariables);
132 });
133 }
134
139 template<std::size_t otherId, class JacRow, class GridVariables,
140 typename std::enable_if_t<(otherId == id), int> = 0>
141 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
142 const LocalResidualValues& res, GridVariables& gridVariables)
143 {}
144
148 template<std::size_t otherId, class JacRow, class GridVariables,
149 typename std::enable_if_t<(otherId != id), int> = 0>
150 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
151 const LocalResidualValues& res, GridVariables& gridVariables)
152 {
153 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
154 }
155
159 void assembleResidual(SubSolutionVector& res)
160 {
161 this->asImp_().bindLocalViews();
162 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
163 res[globalI] = this->evalLocalResidual()[0]; // forward to the internal implementation
164
165 if constexpr (Problem::enableInternalDirichletConstraints())
166 {
167 const auto applyDirichlet = [&] (const auto& scvI,
168 const auto& dirichletValues,
169 const auto eqIdx,
170 const auto pvIdx)
171 {
172 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
173 };
174
175 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
176 }
177 }
178
182 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
183 {
184 // initialize the residual vector for all scvs in this element
185 ElementResidualVector residual(this->fvGeometry().numScv());
186
187 // evaluate the volume terms (storage + source terms)
188 // forward to the local residual specialized for the discretization methods
189 for (auto&& scv : scvs(this->fvGeometry()))
190 {
191 const auto& curVolVars = elemVolVars[scv];
192 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
193 source *= -Extrusion::volume(scv)*curVolVars.extrusionFactor();
194 residual[scv.indexInElement()] = std::move(source);
195 }
196
197 return residual;
198 }
199
203 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
204 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
205
209 LocalResidualValues evalLocalStorageResidual() const
210 {
211 return this->localResidual().evalStorage(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars())[0];
212 }
213
217 LocalResidualValues evalFluxResidual(const Element& neighbor,
218 const SubControlVolumeFace& scvf) const
219 {
220 const auto& elemVolVars = implicit ? this->curElemVolVars() : this->prevElemVolVars();
221 return this->localResidual().evalFlux(problem(), neighbor, this->fvGeometry(), elemVolVars, this->elemFluxVarsCache(), scvf);
222 }
223
228 {
229 // get some references for convenience
230 const auto& element = this->element();
231 const auto& curSol = this->curSol()[domainId];
232 auto&& fvGeometry = this->fvGeometry();
233 auto&& curElemVolVars = this->curElemVolVars();
234 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
235
236 // bind the caches
237 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
238 fvGeometry.bind(element);
239
240 if (implicit)
241 {
244 if (!this->assembler().isStationaryProblem())
245 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
246 }
247 else
248 {
249 auto& prevElemVolVars = this->prevElemVolVars();
250 const auto& prevSol = this->assembler().prevSol()[domainId];
251
253 prevElemVolVars.bind(element, fvGeometry, prevSol);
255 }
256 }
257
259 const Problem& problem() const
260 { return this->assembler().problem(domainId); }
261
263 CouplingManager& couplingManager()
264 { return couplingManager_; }
265
266private:
267 CouplingManager& couplingManager_;
268};
269
281template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
283
290template<std::size_t id, class TypeTag, class Assembler>
291class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
292: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
294{
295 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
296 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
298
301
303 using FVElementGeometry = typename GridGeometry::LocalView;
304 using GridView = typename GridGeometry::GridView;
305 using Element = typename GridView::template Codim<0>::Entity;
306
308 enum { dim = GridView::dimension };
309
310 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
311 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
312 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
313 static constexpr auto domainI = Dune::index_constant<id>();
314
315public:
317
324 template<class JacobianMatrixDiagBlock, class GridVariables>
325 LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
326 {
328 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
329 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
330 // actual element. In the actual element we evaluate the derivative of the entire residual. //
332
333 // get some aliases for convenience
334 const auto& element = this->element();
335 const auto& fvGeometry = this->fvGeometry();
336 const auto& gridGeometry = fvGeometry.gridGeometry();
337 auto&& curElemVolVars = this->curElemVolVars();
338 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
339
340 // get stencil informations
341 const auto globalI = gridGeometry.elementMapper().index(element);
342 const auto& connectivityMap = gridGeometry.connectivityMap();
343 const auto numNeighbors = connectivityMap[globalI].size();
344
345 // container to store the neighboring elements
346 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
347 neighborElements.resize(numNeighbors);
348
349 // assemble the undeflected residual
351 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
352 origResiduals[0] = this->evalLocalResidual()[0];
353
354 // lambda for convenient evaluation of the fluxes across scvfs in the neighbors
355 // if the neighbor is a ghost we don't want to add anything to their residual
356 // so we return 0 and omit computing the flux
357 const auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf)
358 {
359 if (neighbor.partitionType() == Dune::GhostEntity)
360 return LocalResidualValues(0.0);
361 else
362 return this->evalFluxResidual(neighbor, scvf);
363 };
364
365 // get the elements in which we need to evaluate the fluxes
366 // and calculate these in the undeflected state
367 unsigned int j = 1;
368 for (const auto& dataJ : connectivityMap[globalI])
369 {
370 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
371 for (const auto scvfIdx : dataJ.scvfsJ)
372 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
373
374 ++j;
375 }
376
377 // reference to the element's scv (needed later) and corresponding vol vars
378 const auto& scv = fvGeometry.scv(globalI);
379 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
380
381 // save a copy of the original privars and vol vars in order
382 // to restore the original solution after deflection
383 const auto& curSol = this->curSol()[domainI];
384 const auto origPriVars = curSol[globalI];
385 const auto origVolVars = curVolVars;
386
387 // element solution container to be deflected
388 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
389
390 // derivatives in the neighbors with repect to the current elements
391 // in index 0 we save the derivative of the element residual with respect to it's own dofs
392 Residuals partialDerivs(numNeighbors + 1);
393
394 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
395 {
396 partialDerivs = 0.0;
397
398 auto evalResiduals = [&](Scalar priVar)
399 {
400 Residuals partialDerivsTmp(numNeighbors + 1);
401 partialDerivsTmp = 0.0;
402 // update the volume variables and the flux var cache
403 elemSol[0][pvIdx] = priVar;
404 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx);
405 curVolVars.update(elemSol, this->problem(), element, scv);
406 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
407 if (enableGridFluxVarsCache)
408 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
409
410 // calculate the residual with the deflected primary variables
411 partialDerivsTmp[0] = this->evalLocalResidual()[0];
412
413 // calculate the fluxes in the neighbors with the deflected primary variables
414 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
415 for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
416 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
417
418 return partialDerivsTmp;
419 };
420
421 // derive the residuals numerically
422 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
423 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
424 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals,
425 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
426
427 // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the
428 // current primary variable and a 0 elsewhere. As we always solve for a delta of the
429 // solution with repect to the initial one, this results in a delta of 0 for ghosts.
430 if (this->elementIsGhost())
431 {
432 partialDerivs[0] = 0.0;
433 partialDerivs[0][pvIdx] = 1.0;
434 }
435
436 // restore the original state of the scv's volume variables
437 curVolVars = origVolVars;
438
439 // restore the current element solution
440 elemSol[0][pvIdx] = origPriVars[pvIdx];
441
442 // restore the undeflected state of the coupling context
443 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx);
444
445 // add the current partial derivatives to the global jacobian matrix
446 if constexpr (Problem::enableInternalDirichletConstraints())
447 {
448 // check if own element has internal Dirichlet constraint
449 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
450 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
451
452 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
453 {
454 if (internalDirichletConstraintsOwnElement[eqIdx])
455 {
456 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
457 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
458 }
459 else
460 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
461 }
462
463 // off-diagonal entries
464 j = 1;
465 for (const auto& dataJ : connectivityMap[globalI])
466 {
467 const auto& neighborElement = neighborElements[j-1];
468 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
469 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
470
471 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
472 {
473 if (internalDirichletConstraintsNeighbor[eqIdx])
474 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
475 else
476 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
477 }
478
479 ++j;
480 }
481 }
482 else // no internal Dirichlet constraints specified
483 {
484 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
485 {
486 // the diagonal entries
487 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
488
489 // off-diagonal entries
490 j = 1;
491 for (const auto& dataJ : connectivityMap[globalI])
492 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
493 }
494 }
495 }
496
497 // restore original state of the flux vars cache in case of global caching.
498 // This has to be done in order to guarantee that everything is in an undeflected
499 // state before the assembly of another element is called. In the case of local caching
500 // this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
501 // We only have to do this for the last primary variable, for all others the flux var cache
502 // is updated with the correct element volume variables before residual evaluations
503 if (enableGridFluxVarsCache)
504 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
505
506 // compute potential additional derivatives causes by the coupling
507 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
508 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origElementResidual, A, gridVariables);
509
510 // return the original residual
511 return origResiduals[0];
512 }
513
518 template<std::size_t otherId, class JacobianBlock, class GridVariables>
519 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
520 const LocalResidualValues& res, GridVariables& gridVariables)
521 {
523 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
525
526 // get some aliases for convenience
527 const auto& element = this->element();
528 const auto& fvGeometry = this->fvGeometry();
529 const auto& gridGeometry = fvGeometry.gridGeometry();
530 auto&& curElemVolVars = this->curElemVolVars();
531 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
532
533 // get stencil informations
534 const auto globalI = gridGeometry.elementMapper().index(element);
535 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
536 const auto& curSolJ = this->curSol()[domainJ];
537
538 // make sure the flux variables cache does not contain any artifacts from prior differentiation
539 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
540
541 // convenience lambda for call to update self
542 auto updateCoupledVariables = [&] ()
543 {
544 // Update ourself after the context has been modified. Depending on the
545 // type of caching, other objects might have to be updated. All ifs can be optimized away.
546 if (enableGridFluxVarsCache)
547 {
548 if (enableGridVolVarsCache)
549 {
550 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
551 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
552 }
553 else
554 {
555 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
556 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
557 }
558 }
559 else
560 {
561 if (enableGridVolVarsCache)
562 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
563 else
564 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
565 }
566 };
567
568 for (const auto& globalJ : stencil)
569 {
570 // undeflected privars and privars to be deflected
571 const auto origPriVarsJ = curSolJ[globalJ];
572 auto priVarsJ = origPriVarsJ;
573
574 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
575
576 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
577 {
578 auto evalCouplingResidual = [&](Scalar priVar)
579 {
580 // update the volume variables and the flux var cache
581 priVarsJ[pvIdx] = priVar;
582 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
583 updateCoupledVariables();
584 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
585 };
586
587 // derive the residuals numerically
588 LocalResidualValues partialDeriv(0.0);
589 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
590 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
591 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
592 NumericDifferentiation::partialDerivative(evalCouplingResidual, priVarsJ[pvIdx], partialDeriv, origResidual,
593 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
594
595 // add the current partial derivatives to the global jacobian matrix
596 if constexpr (Problem::enableInternalDirichletConstraints())
597 {
598 const auto& scv = fvGeometry.scv(globalI);
599 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
600 if (internalDirichletConstraints.none())
601 {
602 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
603 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
604 }
605 else
606 {
607 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
608 {
609 if (internalDirichletConstraints[eqIdx])
610 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
611 else
612 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
613 }
614 }
615 }
616 else
617 {
618 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
619 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
620 }
621
622 // restore the current element solution
623 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
624
625 // restore the undeflected state of the coupling context
626 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
627 }
628
629 // restore original state of the flux vars cache and/or vol vars in case of global caching.
630 // This has to be done in order to guarantee that the undeflected residual computation done
631 // above is correct when jumping to the next coupled dof of domainJ.
632 updateCoupledVariables();
633 }
634 }
635};
636
643template<std::size_t id, class TypeTag, class Assembler>
644class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
645: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
646 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
647{
648 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
649 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/false>;
650
654
655 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
656 static constexpr auto domainI = Dune::index_constant<id>();
657
658public:
660
671 template<class JacobianMatrixDiagBlock, class GridVariables>
672 LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
673 {
674 if (this->assembler().isStationaryProblem())
675 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
676
677 // assemble the undeflected residual
678 const auto residual = this->evalLocalResidual()[0];
679 const auto storageResidual = this->evalLocalStorageResidual();
680
682 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
683 // neighboring elements all derivatives are zero. For the assembled element only the storage //
684 // derivatives are non-zero. //
686
687 // get some aliases for convenience
688 const auto& element = this->element();
689 const auto& fvGeometry = this->fvGeometry();
690 const auto& gridGeometry = fvGeometry.gridGeometry();
691 auto&& curElemVolVars = this->curElemVolVars();
692
693 // reference to the element's scv (needed later) and corresponding vol vars
694 const auto globalI = gridGeometry.elementMapper().index(element);
695 const auto& scv = fvGeometry.scv(globalI);
696 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
697
698 // save a copy of the original privars and vol vars in order
699 // to restore the original solution after deflection
700 const auto& curSol = this->curSol()[domainI];
701 const auto origPriVars = curSol[globalI];
702 const auto origVolVars = curVolVars;
703
704 // element solution container to be deflected
705 auto elemSol = elementSolution(element, curSol, gridGeometry);
706
707 // derivatives in the neighbors with repect to the current elements
708 LocalResidualValues partialDeriv;
709 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
710 {
711 // reset derivatives of element dof with respect to itself
712 partialDeriv = 0.0;
713
714 auto evalStorage = [&](Scalar priVar)
715 {
716 // update the volume variables and calculate
717 // the residual with the deflected primary variables
718 elemSol[0][pvIdx] = priVar;
719 curVolVars.update(elemSol, this->problem(), element, scv);
720 return this->evalLocalStorageResidual();
721 };
722
723 // for non-ghosts compute the derivative numerically
724 if (!this->elementIsGhost())
725 {
726 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
727 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
728 NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual,
729 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
730 }
731
732 // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else
733 // as we always solve for a delta of the solution with repect to the initial solution this
734 // results in a delta of zero for ghosts
735 else partialDeriv[pvIdx] = 1.0;
736
737 // restore the original state of the scv's volume variables
738 curVolVars = origVolVars;
739
740 // restore the current element solution
741 elemSol[0][pvIdx] = origPriVars[pvIdx];
742
743 // add the current partial derivatives to the global jacobian matrix
744 // only diagonal entries for explicit jacobians
745 if constexpr (Problem::enableInternalDirichletConstraints())
746 {
747 // check if own element has internal Dirichlet constraint
748 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
749 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
750
751 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
752 {
753 if (internalDirichletConstraints[eqIdx])
754 {
755 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
756 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
757 }
758 else
759 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
760 }
761 }
762 else
763 {
764 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
765 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
766 }
767 }
768
769 // return the original residual
770 return residual;
771 }
772
779 template<std::size_t otherId, class JacobianBlock, class GridVariables>
780 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
781 const LocalResidualValues& res, GridVariables& gridVariables)
782 {}
783};
784
790template<std::size_t id, class TypeTag, class Assembler>
791class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
792: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
793 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
794{
795 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>;
796 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
800 using Element = typename GridView::template Codim<0>::Entity;
802
804 enum { dim = GridView::dimension };
805
806 static constexpr auto domainI = Dune::index_constant<id>();
807
808public:
810
817 template<class JacobianMatrixDiagBlock, class GridVariables>
818 LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
819 {
820 // treat ghost separately, we always want zero update for ghosts
821 if (this->elementIsGhost())
822 {
823 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(this->element());
824 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
825 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
826
827 // return zero residual
828 return LocalResidualValues(0.0);
829 }
830
831 // get some aliases for convenience
832 const auto& problem = this->problem();
833 const auto& element = this->element();
834 const auto& fvGeometry = this->fvGeometry();
835 const auto& curElemVolVars = this->curElemVolVars();
836 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
837
838 // get reference to the element's current vol vars
839 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(element);
840 const auto& scv = fvGeometry.scv(globalI);
841 const auto& volVars = curElemVolVars[scv];
842
843 // if the problem is instationary, add derivative of storage term
844 if (!this->assembler().isStationaryProblem())
845 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
846
847 // add source term derivatives
848 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
849
850 // add flux derivatives for each scvf
851 for (const auto& scvf : scvfs(fvGeometry))
852 {
853 // inner faces
854 if (!scvf.boundary())
855 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
856
857 // boundary faces
858 else
859 {
860 const auto& bcTypes = problem.boundaryTypes(element, scvf);
861
862 // add Dirichlet boundary flux derivatives
863 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
864 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
865
866 // add Robin ("solution dependent Neumann") boundary flux derivatives
867 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
868 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
869
870 else
871 DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
872 }
873 }
874
875 if constexpr (Problem::enableInternalDirichletConstraints())
876 {
877 // check if own element has internal Dirichlet constraint
878 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
879 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
880
881 auto residual = this->evalLocalResidual()[0];
882
883 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
884 {
885 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
886 {
887 if (internalDirichletConstraints[eqIdx])
888 {
889 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
890 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
891
892 // inner faces
893 for (const auto& scvf : scvfs(fvGeometry))
894 if (!scvf.boundary())
895 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
896 }
897 }
898 }
899 return residual;
900 }
901 else
902 // return element residual
903 return this->evalLocalResidual()[0];
904 }
905
912 template<std::size_t otherId, class JacobianBlock, class GridVariables>
913 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
914 const LocalResidualValues& res, GridVariables& gridVariables)
915 {
917 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
919
920 // get some aliases for convenience
921 const auto& element = this->element();
922 const auto& fvGeometry = this->fvGeometry();
923 const auto& gridGeometry = fvGeometry.gridGeometry();
924 auto&& curElemVolVars = this->curElemVolVars();
925 // auto&& elemFluxVarsCache = this->elemFluxVarsCache();
926
927 // get stencil informations
928 const auto globalI = gridGeometry.elementMapper().index(element);
929 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
930
931 for (const auto globalJ : stencil)
932 {
933 const auto& elementJ = this->assembler().gridGeometry(domainJ).element(globalJ);
934 this->couplingManager().addCouplingDerivatives(A[globalI][globalJ], domainI, element, fvGeometry, curElemVolVars, domainJ, elementJ);
935
936 // add the current partial derivatives to the global jacobian matrix
937 if constexpr (Problem::enableInternalDirichletConstraints())
938 {
939 const auto& scv = fvGeometry.scv(globalI);
940 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
941 if (internalDirichletConstraints.any())
942 {
943 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
944 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
945 if (internalDirichletConstraints[eqIdx])
946 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
947 }
948 }
949 }
950 }
951};
952
953} // end namespace Dumux
954
955#endif
An adapter class for local assemblers using numeric differentiation.
An enum class to define various differentiation methods available in order to compute the derivatives...
Element solution classes and factory functions.
Helper classes to compute the integration elements.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A helper to deduce a vector with the same size as numbers of equations.
A arithmetic block vector type based on DUNE's reserved vector.
A class for numeric differentiation.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethod::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Definition: common/pdesolver.hh:36
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
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
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with repect to a function parameter.
Definition: numericdifferentiation.hh:61
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:38
Definition: multidomain/couplingmanager.hh:46
A base class for all multidomain local assemblers.
Definition: subdomaincclocalassembler.hh:61
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomaincclocalassembler.hh:227
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: subdomaincclocalassembler.hh:97
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomaincclocalassembler.hh:203
LocalResidualValues evalFluxResidual(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes depending on the chose time discretization scheme.
Definition: subdomaincclocalassembler.hh:217
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomaincclocalassembler.hh:118
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const LocalResidualValues &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomaincclocalassembler.hh:141
const Problem & problem() const
return reference to the underlying problem
Definition: subdomaincclocalassembler.hh:259
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomaincclocalassembler.hh:263
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomaincclocalassembler.hh:90
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
the local residual type of this domain
Definition: subdomaincclocalassembler.hh:92
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomaincclocalassembler.hh:182
LocalResidualValues evalLocalStorageResidual() const
Evaluates the storage terms within the element.
Definition: subdomaincclocalassembler.hh:209
void assembleResidual(SubSolutionVector &res)
Assemble the residual only.
Definition: subdomaincclocalassembler.hh:159
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time disc...
Definition: subdomaincclocalassembler.hh:294
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:519
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:325
Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time disc...
Definition: subdomaincclocalassembler.hh:647
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:672
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the coupling derivatives with respect to the given element and adds them to the global matri...
Definition: subdomaincclocalassembler.hh:780
Cell-centered scheme local assembler using analytic differentiation and implicit time discretization.
Definition: subdomaincclocalassembler.hh:794
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const LocalResidualValues &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:913
LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomaincclocalassembler.hh:818
Declares all properties used in Dumux.