3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fclocalassembler.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 *****************************************************************************/
25#ifndef DUMUX_FC_LOCAL_ASSEMBLER_HH
26#define DUMUX_FC_LOCAL_ASSEMBLER_HH
27
28#include <dune/grid/common/gridenums.hh>
29
39
40namespace Dumux {
41
42namespace Detail {
43
45{
46 template<class... Args>
47 constexpr void operator()(Args&&...) const {}
48};
49
50template<class T, class Default>
51using NonVoidOrDefault_t = std::conditional_t<!std::is_same_v<T, void>, T, Default>;
52
53} // end namespace Detail
54
64template<class TypeTag, class Assembler, class Implementation, bool implicit>
65class FaceCenteredLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
66{
73
74 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
75
76public:
77
78 using ParentType::ParentType;
79
81 {
83 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
84 }
85
90 template <class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::NoOpFunctor>
91 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
92 const PartialReassembler* partialReassembler,
93 const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction{})
94 {
95 static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(),
96 "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!");
97
98 this->asImp_().bindLocalViews();
99 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
100 const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element());
101 if (partialReassembler
102 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
103 {
104 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
105 for (const auto& scv : scvs(this->fvGeometry()))
106 res[scv.dofIndex()] += residual[scv.localDofIndex()];
107
108 // assemble the coupling blocks for coupled models (does nothing if not coupled)
109 maybeAssembleCouplingBlocks(residual);
110 }
111 else if (!this->elementIsGhost())
112 {
113 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
114
115 if (this->element().partitionType() == Dune::InteriorEntity)
116 {
117 for (const auto& scv : scvs(this->fvGeometry()))
118 res[scv.dofIndex()] += residual[scv.localDofIndex()];
119 }
120 else
121 {
122 // handle residual and matrix entries for parallel runs
123 for (const auto& scv : scvs(this->fvGeometry()))
124 {
125 const auto& facet = this->element().template subEntity <1> (scv.indexInElement());
126 // make sure that the residual at border entities is consistent by adding the
127 // the contribution from the neighboring overlap element's scv
128 if (facet.partitionType() == Dune::BorderEntity)
129 res[scv.dofIndex()] += residual[scv.localDofIndex()];
130
131 // set the matrix entries of all DOFs within the overlap region (except the border DOF)
132 // to 1.0 and the residual entries to 0.0
133 else
134 {
135 const auto idx = scv.dofIndex();
136 jac[idx][idx] = 0.0;
137 for (int i = 0; i < jac[idx][idx].size(); ++i)
138 jac[idx][idx][i][i] = 1.0;
139 res[idx] = 0;
140 }
141 }
142 }
143
144 // assemble the coupling blocks for coupled models (does nothing if not coupled)
145 maybeAssembleCouplingBlocks(residual);
146 }
147 else
148 DUNE_THROW(Dune::NotImplemented, "Ghost elements not supported");
149
150
151 auto applyDirichlet = [&] (const auto& scvI,
152 const auto& dirichletValues,
153 const auto eqIdx,
154 const auto pvIdx)
155 {
156 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
157
158 auto& row = jac[scvI.dofIndex()];
159 for (auto col = row.begin(); col != row.end(); ++col)
160 row[col.index()][eqIdx] = 0.0;
161
162 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
163
164 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof TODO periodic
165 // if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
166 // {
167 // const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
168 // res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
169 // const auto end = jac[periodicDof].end();
170 // for (auto it = jac[periodicDof].begin(); it != end; ++it)
171 // (*it) = periodicDof != it.index() ? 0.0 : 1.0;
172 // }
173 };
174
175 this->asImp_().enforceDirichletConstraints(applyDirichlet);
176 }
177
182 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
183 {
184 this->asImp_().bindLocalViews();
185 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
186
187 auto applyDirichlet = [&] (const auto& scvI,
188 const auto& dirichletValues,
189 const auto eqIdx,
190 const auto pvIdx)
191 {
192 auto& row = jac[scvI.dofIndex()];
193 for (auto col = row.begin(); col != row.end(); ++col)
194 row[col.index()][eqIdx] = 0.0;
195
196 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
197 };
198
199 this->asImp_().enforceDirichletConstraints(applyDirichlet);
200 }
201
205 void assembleResidual(SolutionVector& res)
206 {
207 this->asImp_().bindLocalViews();
208 const auto residual = this->evalLocalResidual();
209
210 for (const auto& scv : scvs(this->fvGeometry()))
211 res[scv.dofIndex()] += residual[scv.localDofIndex()];
212
213 auto applyDirichlet = [&] (const auto& scvI,
214 const auto& dirichletValues,
215 const auto eqIdx,
216 const auto pvIdx)
217 {
218 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[eqIdx] - dirichletValues[pvIdx];
219 };
220
221 this->asImp_().enforceDirichletConstraints(applyDirichlet);
222 }
223
227 template<typename ApplyFunction>
228 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
229 {
230 // enforce Dirichlet boundary conditions
231 this->asImp_().evalDirichletBoundaries(applyDirichlet);
232 // take care of internal Dirichlet constraints (if enabled)
233 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
234 }
235
239 template< typename ApplyDirichletFunctionType >
240 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
241 {
242 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
243 // and set the residual to (privar - dirichletvalue)
244 if (this->elemBcTypes().hasDirichlet())
245 {
246 for (const auto& scvf : scvfs(this->fvGeometry()))
247 {
248 if (scvf.isFrontal() && scvf.boundary())
249 {
250 const auto bcTypes = this->elemBcTypes()[scvf.localIndex()];
251 if (bcTypes.hasDirichlet())
252 {
253 const auto& scv = this->fvGeometry().scv(scvf.insideScvIdx());
254 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvf);
255
256 // set the Dirichlet conditions in residual and jacobian
257 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
258 {
259 for (int pvIdx = 0; pvIdx < GridView::dimension; ++pvIdx)
260 {
261 if (bcTypes.isDirichlet(pvIdx) && pvIdx == scv.dofAxis()) // TODO?
262 applyDirichlet(scv, dirichletValues, eqIdx, pvIdx);
263 }
264 }
265 }
266 }
267 }
268 }
269 }
270
275 template<class... Args>
276 void maybeUpdateCouplingContext(Args&&...) {}
277
282 template<class... Args>
284};
285
295template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
297
303template<class TypeTag, class Assembler, class Implementation>
304class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
306 TypeTag, Assembler,
307 Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
308 /*implicit=*/true
309>
310{
314 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
316 using FVElementGeometry = typename GridGeometry::LocalView;
317 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
318 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
322
323 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
325
326public:
327
329 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
331
338 template <class PartialReassembler = DefaultPartialReassembler>
339 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
340 const PartialReassembler* partialReassembler = nullptr)
341 {
342 // get some aliases for convenience
343 const auto& problem = this->asImp_().problem();
344 const auto& element = this->element();
345 const auto& fvGeometry = this->fvGeometry();
346 const auto& curSol = this->asImp_().curSol();
347 auto&& curElemVolVars = this->curElemVolVars();
348
349 // get the vector of the actual element residuals
350 const auto origResiduals = this->evalLocalResidual();
351
353 // //
354 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
355 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
356 // actual element. In the actual element we evaluate the derivative of the entire residual. //
357 // //
359
360 // one residual per element facet
361 const auto numElementResiduals = fvGeometry.numScv();
362
363 // create the vector storing the partial derivatives
364 ElementResidualVector partialDerivs(numElementResiduals);
365
366 const auto evalSource = [&](ElementResidualVector& residual, const SubControlVolume& scv)
367 {
368 this->localResidual().evalSource(residual, problem, element, fvGeometry, curElemVolVars, scv);
369 };
370
371 const auto evalStorage = [&](ElementResidualVector& residual, const SubControlVolume& scv)
372 {
373 this->localResidual().evalStorage(residual, problem, element, fvGeometry, this->prevElemVolVars(), curElemVolVars, scv);
374 };
375
376 const auto evalFlux = [&](ElementResidualVector& residual, const SubControlVolumeFace& scvf)
377 {
378 if (!scvf.processorBoundary())
379 this->localResidual().evalFlux(residual, problem, element, fvGeometry, curElemVolVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf);
380 };
381
382 const auto evalDerivative = [&] (const auto& scvI, const auto& scvJ)
383 {
384 // derivative w.r.t. own DOF
385 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
386 {
387 partialDerivs = 0.0;
388 const auto& otherElement = fvGeometry.gridGeometry().element(scvJ.elementIndex());
389 auto otherElemSol = elementSolution(otherElement, curSol, fvGeometry.gridGeometry()); // TODO allow selective creation of elemsol (for one scv)
390 auto& curOtherVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ);
391 const VolumeVariables origOtherVolVars(curOtherVolVars);
392
393 auto evalResiduals = [&](Scalar priVar)
394 {
395 // update the volume variables and compute element residual
396 otherElemSol[scvJ.localDofIndex()][pvIdx] = priVar;
397 curOtherVolVars.update(otherElemSol, problem, otherElement, scvJ);
398 this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx);
399
400 ElementResidualVector residual(numElementResiduals);
401 residual = 0;
402
403 evalSource(residual, scvI);
404
405 if (!this->assembler().isStationaryProblem())
406 evalStorage(residual, scvI);
407
408 for (const auto& scvf : scvfs(fvGeometry, scvI))
409 evalFlux(residual, scvf);
410
411 return residual;
412 };
413
414 // derive the residuals numerically
415 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
416 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
417 NumericDifferentiation::partialDerivative(evalResiduals, otherElemSol[scvJ.localDofIndex()][pvIdx], partialDerivs, origResiduals,
418 eps_(otherElemSol[scvJ.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
419
420 const auto updateJacobian = [&]()
421 {
422 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
423 {
424 // A[i][col][eqIdx][pvIdx] is the rate of change of
425 // the residual of equation 'eqIdx' at dof 'i'
426 // depending on the primary variable 'pvIdx' at dof
427 // 'col'.
428 A[scvI.dofIndex()][scvJ.dofIndex()][eqIdx][pvIdx] += partialDerivs[scvI.localDofIndex()][eqIdx];
429 }
430 };
431
432 using GeometryHelper = typename std::decay_t<decltype(fvGeometry.gridGeometry())>::GeometryHelper;
433 using LocalIntersectionMapper = typename std::decay_t<decltype(fvGeometry.gridGeometry())>::LocalIntersectionMapper;
434 LocalIntersectionMapper localIsMapper;
435
436 const bool isParallel = fvGeometry.gridGeometry().gridView().comm().size() > 1;
437 if (isParallel)
438 localIsMapper.update(fvGeometry.gridGeometry().gridView(), element);
439
440 if (element.partitionType() == Dune::InteriorEntity)
441 updateJacobian();
442 else
443 {
444 const auto localIdxI = scvI.indexInElement();
445 const auto localIdxJ = scvJ.indexInElement();
446
447 const auto& facetI = GeometryHelper::facet(localIsMapper.refToRealIdx(localIdxI), element);
448 // add contribution of opposite scv lying within the overlap/ghost zone
449 if (facetI.partitionType() == Dune::BorderEntity &&
450 (localIdxJ == GeometryHelper::localOppositeIdx(localIdxI) || scvJ.dofIndex() == scvI.dofIndex()))
451 updateJacobian();
452 }
453
454 if (isParallel && element.partitionType() == Dune::InteriorEntity)
455 {
456 const auto localIdxI = scvI.indexInElement();
457 const auto& facetI = GeometryHelper::facet(localIsMapper.refToRealIdx(localIdxI), element);
458 if (facetI.partitionType() == Dune::BorderEntity)
459 {
460 for (const auto& scvf : scvfs(fvGeometry, scvI))
461 {
462 if (scvf.isFrontal() || scvf.boundary())
463 continue;
464
465 // parallel scvs TODO drawing
466 if (scvf.outsideScvIdx() == scvJ.index())
467 updateJacobian();
468 else
469 {
470 // normal scvs
471 const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(scvf);
472 if (orthogonalScvf.boundary())
473 continue;
474
475 if (orthogonalScvf.insideScvIdx() == scvJ.index() || orthogonalScvf.outsideScvIdx() == scvJ.index())
476 updateJacobian();
477 }
478 }
479 }
480 }
481
482 // restore the original state of the scv's volume variables
483 curOtherVolVars = origOtherVolVars;
484
485 // restore the original element solution
486 otherElemSol[scvJ.localDofIndex()][pvIdx] = curSol[scvJ.dofIndex()][pvIdx];
487 this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx);
488 // TODO additional dof dependencies
489 }
490 };
491
492 // calculation of the derivatives
493 for (const auto& scvI : scvs(fvGeometry))
494 {
495 // derivative w.r.t. own DOFs
496 evalDerivative(scvI, scvI);
497
498 // derivative w.r.t. other DOFs
499 const auto& otherScvIndices = fvGeometry.gridGeometry().connectivityMap()[scvI.index()];
500 for (const auto globalJ : otherScvIndices)
501 evalDerivative(scvI, fvGeometry.scv(globalJ));
502 }
503
504 // evaluate additional derivatives that might arise from the coupling
505 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
506
507 return origResiduals;
508 }
509};
510
511
517template<class TypeTag, class Assembler, class Implementation>
518class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation>
520 TypeTag, Assembler,
521 Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
522 /*implicit=*/false
523>
524{
528 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
532
534
535public:
537 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
539
546 template <class PartialReassembler = DefaultPartialReassembler>
547 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
548 const PartialReassembler* partialReassembler = nullptr)
549 {
550 if (partialReassembler)
551 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
552
553 // get some aliases for convenience
554 const auto& problem = this->asImp_().problem();
555 const auto& element = this->element();
556 const auto& fvGeometry = this->fvGeometry();
557 const auto& curSol = this->asImp_().curSol();
558 auto&& curElemVolVars = this->curElemVolVars();
559
560 // get the vector of the actual element residuals
561 const auto origResiduals = this->evalLocalResidual();
562 const auto origStorageResiduals = this->evalLocalStorageResidual();
563
565 // //
566 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
567 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
568 // actual element. In the actual element we evaluate the derivative of the entire residual. //
569 // //
571
572 // create the element solution
573 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
574
575 // create the vector storing the partial derivatives
576 ElementResidualVector partialDerivs(fvGeometry.numScv());
577
578 // calculation of the derivatives
579 for (const auto& scv : scvs(fvGeometry))
580 {
581 // dof index and corresponding actual pri vars
582 const auto dofIdx = scv.dofIndex();
583 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
584 const VolumeVariables origVolVars(curVolVars);
585
586 // calculate derivatives w.r.t to the privars at the dof at hand
587 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
588 {
589 partialDerivs = 0.0;
590
591 auto evalStorage = [&](Scalar priVar)
592 {
593 // auto partialDerivsTmp = partialDerivs;
594 elemSol[scv.localDofIndex()][pvIdx] = priVar;
595 curVolVars.update(elemSol, problem, element, scv);
596 return this->evalLocalStorageResidual();
597 };
598
599 // derive the residuals numerically
600 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
601 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
602 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
603 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
604
605 // update the global stiffness matrix with the current partial derivatives
606 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
607 {
608 // A[i][col][eqIdx][pvIdx] is the rate of change of
609 // the residual of equation 'eqIdx' at dof 'i'
610 // depending on the primary variable 'pvIdx' at dof
611 // 'col'.
612 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
613 }
614
615 // restore the original state of the scv's volume variables
616 curVolVars = origVolVars;
617
618 // restore the original element solution
619 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
620 }
621 }
622 return origResiduals;
623 }
624};
625
631template<class TypeTag, class Assembler, class Implementation>
632class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation>
634 TypeTag, Assembler,
635 FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>,
636 /*implicit=*/true
637>
638{
644
646
647public:
649 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
651
658 template <class PartialReassembler = DefaultPartialReassembler>
659 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
660 const PartialReassembler* partialReassembler = nullptr)
661 {
662 if (partialReassembler)
663 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
664
665 // get some aliases for convenience
666 const auto& element = this->element();
667 const auto& fvGeometry = this->fvGeometry();
668 const auto& problem = this->asImp_().problem();
669 const auto& curElemVolVars = this->curElemVolVars();
670 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
671
672 // get the vecor of the actual element residuals
673 const auto origResiduals = this->evalLocalResidual();
674
676 // //
677 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
678 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
679 // actual element. In the actual element we evaluate the derivative of the entire residual. //
680 // //
682
683 // calculation of the source and storage derivatives
684 for (const auto& scv : scvs(fvGeometry))
685 {
686 // dof index and corresponding actual pri vars
687 const auto dofIdx = scv.dofIndex();
688 const auto& volVars = curElemVolVars[scv];
689
690 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
691 // only if the problem is instationary we add derivative of storage term
692 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
693 if (!this->assembler().isStationaryProblem())
694 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
695 problem,
696 element,
697 fvGeometry,
698 volVars,
699 scv);
700
701 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
702 // add source term derivatives
703 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
704 problem,
705 element,
706 fvGeometry,
707 volVars,
708 scv);
709 }
710
711 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
712 for (const auto& scvf : scvfs(fvGeometry))
713 {
714 if (!scvf.boundary())
715 {
716 // add flux term derivatives
717 this->localResidual().addFluxDerivatives(A,
718 problem,
719 element,
720 fvGeometry,
721 curElemVolVars,
722 elemFluxVarsCache,
723 scvf);
724 }
725
726 // the boundary gets special treatment to simplify
727 // for the user
728 else
729 {
730 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
731 if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
732 {
733 // add flux term derivatives
734 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
735 problem,
736 element,
737 fvGeometry,
738 curElemVolVars,
739 elemFluxVarsCache,
740 scvf);
741 }
742 }
743 }
744
745 return origResiduals;
746 }
747};
748
754template<class TypeTag, class Assembler, class Implementation>
755class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation>
757 TypeTag, Assembler,
758 FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>,
759 /*implicit=*/false
760>
761{
767
769
770public:
772 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
774
781 template <class PartialReassembler = DefaultPartialReassembler>
782 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
783 const PartialReassembler* partialReassembler = nullptr)
784 {
785 if (partialReassembler)
786 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
787
788 // get some aliases for convenience
789 const auto& element = this->element();
790 const auto& fvGeometry = this->fvGeometry();
791 const auto& problem = this->asImp_().problem();
792 const auto& curElemVolVars = this->curElemVolVars();
793
794 // get the vector of the actual element residuals
795 const auto origResiduals = this->evalLocalResidual();
796
798 // //
799 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
800 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
801 // actual element. In the actual element we evaluate the derivative of the entire residual. //
802 // //
804
805 // calculation of the source and storage derivatives
806 for (const auto& scv : scvs(fvGeometry))
807 {
808 // dof index and corresponding actual pri vars
809 const auto dofIdx = scv.dofIndex();
810 const auto& volVars = curElemVolVars[scv];
811
812 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
813 // only if the problem is instationary we add derivative of storage term
814 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
815 problem,
816 element,
817 fvGeometry,
818 volVars,
819 scv);
820 }
821
822 return origResiduals;
823 }
824};
825
826} // end namespace Dumux
827
828#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
Detects which entries in the Jacobian have to be recomputed.
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::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
@ green
does not need to be reassembled
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
std::conditional_t<!std::is_same_v< T, void >, T, Default > NonVoidOrDefault_t
Definition: fclocalassembler.hh:51
Definition: fclocalassembler.hh:45
constexpr void operator()(Args &&...) const
Definition: fclocalassembler.hh:47
A base class for all local cell-centered assemblers.
Definition: fclocalassembler.hh:66
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: fclocalassembler.hh:228
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: fclocalassembler.hh:205
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fclocalassembler.hh:182
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler, const CouplingFunction &maybeAssembleCouplingBlocks=CouplingFunction{})
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: fclocalassembler.hh:91
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: fclocalassembler.hh:276
void bindLocalViews()
Definition: fclocalassembler.hh:80
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: fclocalassembler.hh:240
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: fclocalassembler.hh:283
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:296
Face-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: fclocalassembler.hh:310
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fclocalassembler.hh:329
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fclocalassembler.hh:339
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fclocalassembler.hh:328
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fclocalassembler.hh:537
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fclocalassembler.hh:536
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fclocalassembler.hh:547
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fclocalassembler.hh:649
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fclocalassembler.hh:659
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fclocalassembler.hh:648
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fclocalassembler.hh:771
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fclocalassembler.hh:772
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fclocalassembler.hh:782
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: fvlocalassemblerbase.hh:185
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
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
const Problem & problem() const
The problem.
Definition: fvlocalassemblerbase.hh:241
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:503
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
Declares all properties used in Dumux.