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