version 3.8
multidomain/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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_CC_LOCAL_ASSEMBLER_HH
16
17#include <dune/common/indices.hh>
18#include <dune/common/reservedvector.hh>
19#include <dune/common/hybridutilities.hh>
20#include <dune/grid/common/gridenums.hh> // for GhostEntity
21#include <dune/istl/matrixindexset.hh>
22
33
34namespace Dumux {
35
47template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool implicit>
48class SubDomainCCLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler,Implementation, implicit>
49{
51
55 using SolutionVector = typename Assembler::SolutionVector;
57
59 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
60 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
61 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
62 using Scalar = typename GridVariables::Scalar;
63
64 using GridGeometry = typename GridVariables::GridGeometry;
65 using FVElementGeometry = typename GridGeometry::LocalView;
66 using SubControlVolume = typename GridGeometry::SubControlVolume;
67 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
68 using Extrusion = Extrusion_t<GridGeometry>;
69 using GridView = typename GridGeometry::GridView;
70 using Element = typename GridView::template Codim<0>::Entity;
71
72 using CouplingManager = typename Assembler::CouplingManager;
73 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
74
75public:
77 static constexpr auto domainId = typename Dune::index_constant<id>();
81 using ParentType::ParentType;
82
84 explicit SubDomainCCLocalAssemblerBase(const Assembler& assembler,
85 const Element& element,
86 const SolutionVector& curSol,
87 CouplingManager& couplingManager)
89 element,
90 curSol,
91 localView(assembler.gridGeometry(domainId)),
92 localView(assembler.gridVariables(domainId).curGridVolVars()),
93 localView(assembler.gridVariables(domainId).prevGridVolVars()),
94 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
96 (element.partitionType() == Dune::GhostEntity))
97 , couplingManager_(couplingManager)
98 {}
99
104 template<class JacobianMatrixRow, class SubResidualVector, class GridVariablesTuple>
105 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubResidualVector& res, GridVariablesTuple& gridVariables)
106 {
107 this->asImp_().bindLocalViews();
108
109 // for the diagonal jacobian block
110 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
111 res[globalI] = this->asImp_().assembleJacobianAndResidualImplInverse(jacRow[domainId], *std::get<domainId>(gridVariables));
112
113 // for the coupling blocks
114 using namespace Dune::Hybrid;
115 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
116 {
117 if (i != id)
118 this->assembleJacobianCoupling(i, jacRow, res[globalI], gridVariables);
119 });
120 }
121
126 template<std::size_t otherId, class JacRow, class GridVariables,
127 typename std::enable_if_t<(otherId == id), int> = 0>
128 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
129 const LocalResidualValues& res, GridVariables& gridVariables)
130 {}
131
135 template<std::size_t otherId, class JacRow, class GridVariables,
136 typename std::enable_if_t<(otherId != id), int> = 0>
137 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
138 const LocalResidualValues& res, GridVariables& gridVariables)
139 {
140 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
141 }
142
146 template<class SubResidualVector>
147 void assembleResidual(SubResidualVector& res)
148 {
149 this->asImp_().bindLocalViews();
150 const auto globalI = this->fvGeometry().gridGeometry().elementMapper().index(this->element());
151 res[globalI] = this->evalLocalResidual()[0]; // forward to the internal implementation
152
153 if constexpr (Problem::enableInternalDirichletConstraints())
154 {
155 const auto applyDirichlet = [&] (const auto& scvI,
156 const auto& dirichletValues,
157 const auto eqIdx,
158 const auto pvIdx)
159 {
160 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
161 };
162
163 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
164 }
165 }
166
170 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
171 {
172 // initialize the residual vector for all scvs in this element
173 ElementResidualVector residual(this->fvGeometry().numScv());
174
175 // evaluate the volume terms (storage + source terms)
176 // forward to the local residual specialized for the discretization methods
177 for (auto&& scv : scvs(this->fvGeometry()))
178 {
179 const auto& curVolVars = elemVolVars[scv];
180 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
181 source *= -Extrusion::volume(this->fvGeometry(), scv)*curVolVars.extrusionFactor();
182 residual[scv.indexInElement()] = std::move(source);
183 }
184
185 return residual;
186 }
187
191 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
192 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
193
197 LocalResidualValues evalLocalStorageResidual() const
198 {
199 return this->localResidual().evalStorage(this->element(), this->fvGeometry(), this->prevElemVolVars(), this->curElemVolVars())[0];
200 }
201
205 LocalResidualValues evalFluxResidual(const Element& neighbor,
206 const SubControlVolumeFace& scvf) const
207 {
208 const auto& elemVolVars = implicit ? this->curElemVolVars() : this->prevElemVolVars();
209 return this->localResidual().evalFlux(problem(), neighbor, this->fvGeometry(), elemVolVars, this->elemFluxVarsCache(), scvf);
210 }
211
216 {
217 // get some references for convenience
218 const auto& element = this->element();
219 const auto& curSol = this->curSol()[domainId];
220 auto&& fvGeometry = this->fvGeometry();
221 auto&& curElemVolVars = this->curElemVolVars();
222 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
223
224 // bind the caches
225 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
226 fvGeometry.bind(element);
227
228 if (implicit)
229 {
232 if (!this->assembler().isStationaryProblem())
233 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
234 }
235 else
236 {
237 auto& prevElemVolVars = this->prevElemVolVars();
238 const auto& prevSol = this->assembler().prevSol()[domainId];
239
241 prevElemVolVars.bind(element, fvGeometry, prevSol);
243 }
244 }
245
247 const Problem& problem() const
248 { return this->assembler().problem(domainId); }
249
251 CouplingManager& couplingManager()
252 { return couplingManager_; }
253
254private:
255 CouplingManager& couplingManager_;
256};
257
269template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
271
278template<std::size_t id, class TypeTag, class Assembler>
279class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
280: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
282{
283 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
284 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
286
289
291 using FVElementGeometry = typename GridGeometry::LocalView;
292 using GridView = typename GridGeometry::GridView;
293 using Element = typename GridView::template Codim<0>::Entity;
294
296 enum { dim = GridView::dimension };
297
300 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
301 static constexpr auto domainI = Dune::index_constant<id>();
302
303public:
305
312 template<class JacobianMatrixDiagBlock, class GridVariables>
313 LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
314 {
316 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
317 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
318 // actual element. In the actual element we evaluate the derivative of the entire residual. //
320
321 // get some aliases for convenience
322 const auto& element = this->element();
323 const auto& fvGeometry = this->fvGeometry();
324 const auto& gridGeometry = fvGeometry.gridGeometry();
325 auto&& curElemVolVars = this->curElemVolVars();
326 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
327
328 // get stencil information
329 const auto globalI = gridGeometry.elementMapper().index(element);
330 const auto& connectivityMap = gridGeometry.connectivityMap();
331 const auto numNeighbors = connectivityMap[globalI].size();
332
333 // container to store the neighboring elements
334 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
335 neighborElements.resize(numNeighbors);
336
337 // assemble the undeflected residual
339 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
340 origResiduals[0] = this->evalLocalResidual()[0];
341
342 // lambda for convenient evaluation of the fluxes across scvfs in the neighbors
343 // if the neighbor is a ghost we don't want to add anything to their residual
344 // so we return 0 and omit computing the flux
345 const auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf)
346 {
347 if (neighbor.partitionType() == Dune::GhostEntity)
348 return LocalResidualValues(0.0);
349 else
350 return this->evalFluxResidual(neighbor, scvf);
351 };
352
353 // get the elements in which we need to evaluate the fluxes
354 // and calculate these in the undeflected state
355 unsigned int j = 1;
356 for (const auto& dataJ : connectivityMap[globalI])
357 {
358 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
359 for (const auto scvfIdx : dataJ.scvfsJ)
360 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
361
362 ++j;
363 }
364
365 // reference to the element's scv (needed later) and corresponding vol vars
366 const auto& scv = fvGeometry.scv(globalI);
367 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
368
369 // save a copy of the original privars and vol vars in order
370 // to restore the original solution after deflection
371 const auto& curSol = this->curSol()[domainI];
372 const auto origPriVars = curSol[globalI];
373 const auto origVolVars = curVolVars;
374
375 // element solution container to be deflected
376 auto elemSol = elementSolution<FVElementGeometry>(origPriVars);
377
378 // derivatives in the neighbors with respect to the current elements
379 // in index 0 we save the derivative of the element residual with respect to it's own dofs
380 Residuals partialDerivs(numNeighbors + 1);
381
382 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
383 {
384 partialDerivs = 0.0;
385
386 auto evalResiduals = [&](Scalar priVar)
387 {
388 Residuals partialDerivsTmp(numNeighbors + 1);
389 partialDerivsTmp = 0.0;
390 // update the volume variables and the flux var cache
391 elemSol[0][pvIdx] = priVar;
392 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx);
393 curVolVars.update(elemSol, this->problem(), element, scv);
394 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
395 if (enableGridFluxVarsCache)
396 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
397
398 // calculate the residual with the deflected primary variables
399 partialDerivsTmp[0] = this->evalLocalResidual()[0];
400
401 // calculate the fluxes in the neighbors with the deflected primary variables
402 for (std::size_t k = 0; k < connectivityMap[globalI].size(); ++k)
403 for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
404 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
405
406 return partialDerivsTmp;
407 };
408
409 // derive the residuals numerically
410 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
411 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
412 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals,
413 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
414
415 // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the
416 // current primary variable and a 0 elsewhere. As we always solve for a delta of the
417 // solution with respect to the initial one, this results in a delta of 0 for ghosts.
418 if (this->elementIsGhost())
419 {
420 partialDerivs[0] = 0.0;
421 partialDerivs[0][pvIdx] = 1.0;
422 }
423
424 // restore the original state of the scv's volume variables
425 curVolVars = origVolVars;
426
427 // restore the current element solution
428 elemSol[0][pvIdx] = origPriVars[pvIdx];
429
430 // restore the undeflected state of the coupling context
431 this->couplingManager().updateCouplingContext(domainI, *this, domainI, globalI, elemSol[0], pvIdx);
432
433 // add the current partial derivatives to the global jacobian matrix
434 if constexpr (Problem::enableInternalDirichletConstraints())
435 {
436 // check if own element has internal Dirichlet constraint
437 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
438 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
439
440 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
441 {
442 if (internalDirichletConstraintsOwnElement[eqIdx])
443 {
444 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
445 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
446 }
447 else
448 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
449 }
450
451 // off-diagonal entries
452 j = 1;
453 for (const auto& dataJ : connectivityMap[globalI])
454 {
455 const auto& neighborElement = neighborElements[j-1];
456 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
457 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
458
459 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
460 {
461 if (internalDirichletConstraintsNeighbor[eqIdx])
462 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
463 else
464 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
465 }
466
467 ++j;
468 }
469 }
470 else // no internal Dirichlet constraints specified
471 {
472 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
473 {
474 // the diagonal entries
475 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
476
477 // off-diagonal entries
478 j = 1;
479 for (const auto& dataJ : connectivityMap[globalI])
480 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
481 }
482 }
483 }
484
485 // restore original state of the flux vars cache in case of global caching.
486 // This has to be done in order to guarantee that everything is in an undeflected
487 // state before the assembly of another element is called. In the case of local caching
488 // this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
489 // We only have to do this for the last primary variable, for all others the flux var cache
490 // is updated with the correct element volume variables before residual evaluations
491 if (enableGridFluxVarsCache)
492 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
493
494 // compute potential additional derivatives causes by the coupling
495 typename ParentType::LocalResidual::ElementResidualVector origElementResidual({origResiduals[0]});
496 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origElementResidual, A, gridVariables);
497
498 // return the original residual
499 return origResiduals[0];
500 }
501
506 template<std::size_t otherId, class JacobianBlock, class GridVariables>
507 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
508 const LocalResidualValues& res, GridVariables& gridVariables)
509 {
511 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
513
514 // get some aliases for convenience
515 const auto& element = this->element();
516 const auto& fvGeometry = this->fvGeometry();
517 const auto& gridGeometry = fvGeometry.gridGeometry();
518 auto&& curElemVolVars = this->curElemVolVars();
519 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
520
521 // get stencil information
522 const auto globalI = gridGeometry.elementMapper().index(element);
523 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
524 const auto& curSolJ = this->curSol()[domainJ];
525
526 // make sure the flux variables cache does not contain any artifacts from prior differentiation
527 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
528
529 // convenience lambda for call to update self
530 auto updateCoupledVariables = [&] ()
531 {
532 // Update ourself after the context has been modified. Depending on the
533 // type of caching, other objects might have to be updated. All ifs can be optimized away.
534 if (enableGridFluxVarsCache)
535 {
536 if (enableGridVolVarsCache)
537 {
538 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
539 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
540 }
541 else
542 {
543 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
544 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
545 }
546 }
547 else
548 {
549 if (enableGridVolVarsCache)
550 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
551 else
552 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
553 }
554 };
555
556 for (const auto& globalJ : stencil)
557 {
558 // undeflected privars and privars to be deflected
559 const auto origPriVarsJ = curSolJ[globalJ];
560 auto priVarsJ = origPriVarsJ;
561
562 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
563
564 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
565 {
566 auto evalCouplingResidual = [&](Scalar priVar)
567 {
568 // update the volume variables and the flux var cache
569 priVarsJ[pvIdx] = priVar;
570 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
571 updateCoupledVariables();
572 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ)[0];
573 };
574
575 // derive the residuals numerically
576 LocalResidualValues partialDeriv(0.0);
577 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
578 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
579 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
580 NumericDifferentiation::partialDerivative(evalCouplingResidual, priVarsJ[pvIdx], partialDeriv, origResidual,
581 epsCoupl(priVarsJ[pvIdx], pvIdx), numDiffMethod);
582
583 // add the current partial derivatives to the global jacobian matrix
584 if constexpr (Problem::enableInternalDirichletConstraints())
585 {
586 const auto& scv = fvGeometry.scv(globalI);
587 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
588 if (internalDirichletConstraints.none())
589 {
590 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
591 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
592 }
593 else
594 {
595 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
596 {
597 if (internalDirichletConstraints[eqIdx])
598 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
599 else
600 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
601 }
602 }
603 }
604 else
605 {
606 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
607 A[globalI][globalJ][eqIdx][pvIdx] += partialDeriv[eqIdx];
608 }
609
610 // restore the current element solution
611 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
612
613 // restore the undeflected state of the coupling context
614 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
615 }
616
617 // restore original state of the flux vars cache and/or vol vars in case of global caching.
618 // This has to be done in order to guarantee that the undeflected residual computation done
619 // above is correct when jumping to the next coupled dof of domainJ.
620 updateCoupledVariables();
621 }
622 }
623};
624
631template<std::size_t id, class TypeTag, class Assembler>
632class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
633: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
634 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false>
635{
636 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
637 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/false>;
638
642
643 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
644 static constexpr auto domainI = Dune::index_constant<id>();
645
646public:
648
659 template<class JacobianMatrixDiagBlock, class GridVariables>
660 LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
661 {
662 if (this->assembler().isStationaryProblem())
663 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
664
665 // assemble the undeflected residual
666 const auto residual = this->evalLocalResidual()[0];
667 const auto storageResidual = this->evalLocalStorageResidual();
668
670 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
671 // neighboring elements all derivatives are zero. For the assembled element only the storage //
672 // derivatives are non-zero. //
674
675 // get some aliases for convenience
676 const auto& element = this->element();
677 const auto& fvGeometry = this->fvGeometry();
678 const auto& gridGeometry = fvGeometry.gridGeometry();
679 auto&& curElemVolVars = this->curElemVolVars();
680
681 // reference to the element's scv (needed later) and corresponding vol vars
682 const auto globalI = gridGeometry.elementMapper().index(element);
683 const auto& scv = fvGeometry.scv(globalI);
684 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
685
686 // save a copy of the original privars and vol vars in order
687 // to restore the original solution after deflection
688 const auto& curSol = this->curSol()[domainI];
689 const auto origPriVars = curSol[globalI];
690 const auto origVolVars = curVolVars;
691
692 // element solution container to be deflected
693 auto elemSol = elementSolution(element, curSol, gridGeometry);
694
695 // derivatives in the neighbors with respect to the current elements
696 LocalResidualValues partialDeriv;
697 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
698 {
699 // reset derivatives of element dof with respect to itself
700 partialDeriv = 0.0;
701
702 auto evalStorage = [&](Scalar priVar)
703 {
704 // update the volume variables and calculate
705 // the residual with the deflected primary variables
706 elemSol[0][pvIdx] = priVar;
707 curVolVars.update(elemSol, this->problem(), element, scv);
708 return this->evalLocalStorageResidual();
709 };
710
711 // for non-ghosts compute the derivative numerically
712 if (!this->elementIsGhost())
713 {
714 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
715 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
716 NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual,
717 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
718 }
719
720 // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else
721 // as we always solve for a delta of the solution with respect to the initial solution this
722 // results in a delta of zero for ghosts
723 else partialDeriv[pvIdx] = 1.0;
724
725 // restore the original state of the scv's volume variables
726 curVolVars = origVolVars;
727
728 // restore the current element solution
729 elemSol[0][pvIdx] = origPriVars[pvIdx];
730
731 // add the current partial derivatives to the global jacobian matrix
732 // only diagonal entries for explicit jacobians
733 if constexpr (Problem::enableInternalDirichletConstraints())
734 {
735 // check if own element has internal Dirichlet constraint
736 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
737 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
738
739 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
740 {
741 if (internalDirichletConstraints[eqIdx])
742 {
743 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
744 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
745 }
746 else
747 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
748 }
749 }
750 else
751 {
752 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
753 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
754 }
755 }
756
757 // return the original residual
758 return residual;
759 }
760
767 template<std::size_t otherId, class JacobianBlock, class GridVariables>
768 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
769 const LocalResidualValues& res, GridVariables& gridVariables)
770 {}
771};
772
779template<std::size_t id, class TypeTag, class Assembler>
780class SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
781: public SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler,
782 SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, true>, true>
783{
784 using ThisType = SubDomainCCLocalAssembler<id, TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>;
785 using ParentType = SubDomainCCLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
789 using Element = typename GridView::template Codim<0>::Entity;
791
793 enum { dim = GridView::dimension };
794
795 static constexpr auto domainI = Dune::index_constant<id>();
796
797public:
799
806 template<class JacobianMatrixDiagBlock, class GridVariables>
807 LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
808 {
809 // treat ghost separately, we always want zero update for ghosts
810 if (this->elementIsGhost())
811 {
812 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(this->element());
813 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
814 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
815
816 // return zero residual
817 return LocalResidualValues(0.0);
818 }
819
820 // get some aliases for convenience
821 const auto& problem = this->problem();
822 const auto& element = this->element();
823 const auto& fvGeometry = this->fvGeometry();
824 const auto& curElemVolVars = this->curElemVolVars();
825 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
826
827 // get reference to the element's current vol vars
828 const auto globalI = this->assembler().gridGeometry(domainI).elementMapper().index(element);
829 const auto& scv = fvGeometry.scv(globalI);
830 const auto& volVars = curElemVolVars[scv];
831
832 // if the problem is instationary, add derivative of storage term
833 if (!this->assembler().isStationaryProblem())
834 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
835
836 // add source term derivatives
837 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
838
839 // add flux derivatives for each scvf
840 for (const auto& scvf : scvfs(fvGeometry))
841 {
842 // inner faces
843 if (!scvf.boundary())
844 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
845
846 // boundary faces
847 else
848 {
849 const auto& bcTypes = problem.boundaryTypes(element, scvf);
850
851 // add Dirichlet boundary flux derivatives
852 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
853 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
854
855 // add Robin ("solution dependent Neumann") boundary flux derivatives
856 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
857 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
858
859 else
860 DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
861 }
862 }
863
864 if constexpr (Problem::enableInternalDirichletConstraints())
865 {
866 // check if own element has internal Dirichlet constraint
867 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
868 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
869
870 auto residual = this->evalLocalResidual()[0];
871
872 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
873 {
874 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
875 {
876 if (internalDirichletConstraints[eqIdx])
877 {
878 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
879 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
880
881 // inner faces
882 for (const auto& scvf : scvfs(fvGeometry))
883 if (!scvf.boundary())
884 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
885 }
886 }
887 }
888 return residual;
889 }
890 else
891 // return element residual
892 return this->evalLocalResidual()[0];
893 }
894
899 template<std::size_t otherId, class JacobianBlock, class GridVariables>
900 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
901 const LocalResidualValues& res, GridVariables& gridVariables)
902 {
904 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
906
907 // get some aliases for convenience
908 const auto& element = this->element();
909 const auto& fvGeometry = this->fvGeometry();
910 const auto& gridGeometry = fvGeometry.gridGeometry();
911 auto&& curElemVolVars = this->curElemVolVars();
912 // auto&& elemFluxVarsCache = this->elemFluxVarsCache();
913
914 // get stencil information
915 const auto globalI = gridGeometry.elementMapper().index(element);
916 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
917
918 for (const auto globalJ : stencil)
919 {
920 const auto& elementJ = this->assembler().gridGeometry(domainJ).element(globalJ);
921 this->couplingManager().addCouplingDerivatives(A[globalI][globalJ], domainI, element, fvGeometry, curElemVolVars, domainJ, elementJ);
922
923 // add the current partial derivatives to the global jacobian matrix
924 if constexpr (Problem::enableInternalDirichletConstraints())
925 {
926 const auto& scv = fvGeometry.scv(globalI);
927 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
928 if (internalDirichletConstraints.any())
929 {
930 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
931 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
932 if (internalDirichletConstraints[eqIdx])
933 A[globalI][globalJ][eqIdx][pvIdx] = 0.0;
934 }
935 }
936 }
937 }
938};
939
940} // end namespace Dumux
941
942#endif
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:48
A base class for all local assemblers.
Definition: assembly/fvlocalassemblerbase.hh:36
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:253
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/fvlocalassemblerbase.hh:56
Implementation & asImp_()
Definition: assembly/fvlocalassemblerbase.hh:297
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: assembly/fvlocalassemblerbase.hh:257
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: assembly/fvlocalassemblerbase.hh:108
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: assembly/fvlocalassemblerbase.hh:249
const Assembler & assembler() const
The assembler.
Definition: assembly/fvlocalassemblerbase.hh:233
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: assembly/fvlocalassemblerbase.hh:261
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition: assembly/fvlocalassemblerbase.hh:55
LocalResidual & localResidual()
The local residual for the current element.
Definition: assembly/fvlocalassemblerbase.hh:265
const Element & element() const
The current element.
Definition: assembly/fvlocalassemblerbase.hh:237
const SolutionVector & curSol() const
The current solution.
Definition: assembly/fvlocalassemblerbase.hh:245
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:49
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:29
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time disc...
Definition: multidomain/subdomaincclocalassembler.hh:282
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: multidomain/subdomaincclocalassembler.hh:507
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: multidomain/subdomaincclocalassembler.hh:313
Cell-centered scheme local assembler using analytic differentiation and implicit time discretization.
Definition: multidomain/subdomaincclocalassembler.hh:783
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: multidomain/subdomaincclocalassembler.hh:900
LocalResidualValues assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: multidomain/subdomaincclocalassembler.hh:807
Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time disc...
Definition: multidomain/subdomaincclocalassembler.hh:635
LocalResidualValues assembleJacobianAndResidualImplInverse(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: multidomain/subdomaincclocalassembler.hh:660
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: multidomain/subdomaincclocalassembler.hh:768
A base class for all multidomain local assemblers.
Definition: multidomain/subdomaincclocalassembler.hh:49
typename ParentType::LocalResidual LocalResidual
the local residual type of this domain
Definition: multidomain/subdomaincclocalassembler.hh:79
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: multidomain/subdomaincclocalassembler.hh:215
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubResidualVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: multidomain/subdomaincclocalassembler.hh:105
SubDomainCCLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
the constructor
Definition: multidomain/subdomaincclocalassembler.hh:84
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: multidomain/subdomaincclocalassembler.hh:191
LocalResidualValues evalFluxResidual(const Element &neighbor, const SubControlVolumeFace &scvf) const
Evaluates the fluxes depending on the chose time discretization scheme.
Definition: multidomain/subdomaincclocalassembler.hh:205
void assembleResidual(SubResidualVector &res)
Assemble the residual only.
Definition: multidomain/subdomaincclocalassembler.hh:147
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: multidomain/subdomaincclocalassembler.hh:128
const Problem & problem() const
return reference to the underlying problem
Definition: multidomain/subdomaincclocalassembler.hh:247
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: multidomain/subdomaincclocalassembler.hh:251
static constexpr auto domainId
export the domain id of this sub-domain
Definition: multidomain/subdomaincclocalassembler.hh:77
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: multidomain/subdomaincclocalassembler.hh:170
LocalResidualValues evalLocalStorageResidual() const
Evaluates the storage terms within the element.
Definition: multidomain/subdomaincclocalassembler.hh:197
The cell-centered scheme multidomain local assembler.
Definition: multidomain/subdomaincclocalassembler.hh:270
Defines all properties used in Dumux.
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.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:25
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:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Definition: common/pdesolver.hh:24
A helper to deduce a vector with the same size as numbers of equations.
A class for numeric differentiation.
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A arithmetic block vector type based on DUNE's reserved vector.