version 3.10.0
Loading...
Searching...
No Matches
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-FileCopyrightText: 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{
61
62 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
63
64public:
65
66 using ParentType::ParentType;
67
69 {
71 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
72 }
73
78 template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::NoOpFunctor>
79 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
80 const PartialReassembler* partialReassembler,
81 const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction{})
82 {
83 static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(),
84 "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!");
85
86 this->asImp_().bindLocalViews();
87 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
88 const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element());
89 if (partialReassembler
90 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
91 {
92 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
93 for (const auto& scv : scvs(this->fvGeometry()))
94 res[scv.dofIndex()] += residual[scv.localDofIndex()];
95
96 // assemble the coupling blocks for coupled models (does nothing if not coupled)
97 maybeAssembleCouplingBlocks(residual);
98 }
99 else if (!this->elementIsGhost())
100 {
101 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
102
103 if (this->element().partitionType() == Dune::InteriorEntity)
104 {
105 for (const auto& scv : scvs(this->fvGeometry()))
106 res[scv.dofIndex()] += residual[scv.localDofIndex()];
107 }
108 else
109 {
110 // handle residual and matrix entries for parallel runs
111 for (const auto& scv : scvs(this->fvGeometry()))
112 {
113 const auto& facet = this->element().template subEntity <1> (scv.indexInElement());
114 // make sure that the residual at border entities is consistent by adding the
115 // the contribution from the neighboring overlap element's scv
116 if (facet.partitionType() == Dune::BorderEntity)
117 res[scv.dofIndex()] += residual[scv.localDofIndex()];
118
119 // set the matrix entries of all DOFs within the overlap region (except the border DOF)
120 // to 1.0 and the residual entries to 0.0
121 else
122 {
123 const auto idx = scv.dofIndex();
124 jac[idx][idx] = 0.0;
125 for (int i = 0; i < jac[idx][idx].size(); ++i)
126 jac[idx][idx][i][i] = 1.0;
127 res[idx] = 0;
128 }
129 }
130 }
131
132 // assemble the coupling blocks for coupled models (does nothing if not coupled)
133 maybeAssembleCouplingBlocks(residual);
134 }
135 else
136 DUNE_THROW(Dune::NotImplemented, "Ghost elements not supported");
137
138
139 auto applyDirichlet = [&] (const auto& scvI,
140 const auto& dirichletValues,
141 const auto eqIdx,
142 const auto pvIdx)
143 {
144 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
145
146 auto& row = jac[scvI.dofIndex()];
147 for (auto col = row.begin(); col != row.end(); ++col)
148 row[col.index()][eqIdx] = 0.0;
149
150 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
151
152 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
153 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
154 {
155 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
156 res[periodicDof][eqIdx] = this->asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx];
157
158 auto& rowP = jac[periodicDof];
159 for (auto col = rowP.begin(); col != rowP.end(); ++col)
160 row[col.index()][eqIdx] = 0.0;
161
162 rowP[periodicDof][eqIdx][pvIdx] = 1.0;
163 }
164 };
165
166 this->asImp_().enforceDirichletConstraints(applyDirichlet);
167 }
168
173 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
174 {
175 this->asImp_().bindLocalViews();
176 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
177
178 auto applyDirichlet = [&] (const auto& scvI,
179 const auto& dirichletValues,
180 const auto eqIdx,
181 const auto pvIdx)
182 {
183 auto& row = jac[scvI.dofIndex()];
184 for (auto col = row.begin(); col != row.end(); ++col)
185 row[col.index()][eqIdx] = 0.0;
186
187 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
188 };
189
190 this->asImp_().enforceDirichletConstraints(applyDirichlet);
191 }
192
196 template<class ResidualVector>
197 void assembleResidual(ResidualVector& res)
198 {
199 this->asImp_().bindLocalViews();
200 const auto residual = this->evalLocalResidual();
201
202 for (const auto& scv : scvs(this->fvGeometry()))
203 res[scv.dofIndex()] += residual[scv.localDofIndex()];
204
205 auto applyDirichlet = [&] (const auto& scvI,
206 const auto& dirichletValues,
207 const auto eqIdx,
208 const auto pvIdx)
209 {
210 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
211 };
212
213 this->asImp_().enforceDirichletConstraints(applyDirichlet);
214 }
215
219 template<typename ApplyFunction>
220 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
221 {
222 // enforce Dirichlet boundary conditions
223 this->asImp_().evalDirichletBoundaries(applyDirichlet);
224 // take care of internal Dirichlet constraints (if enabled)
225 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
226 }
227
231 template< typename ApplyDirichletFunctionType >
232 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
233 {
234 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
235 // and set the residual to (privar - dirichletvalue)
236 if (this->elemBcTypes().hasDirichlet())
237 {
238 for (const auto& scvf : scvfs(this->fvGeometry()))
239 {
240 if (scvf.isFrontal() && scvf.boundary())
241 {
242 const auto bcTypes = this->elemBcTypes()[scvf.localIndex()];
243 if (bcTypes.hasDirichlet())
244 {
245 const auto& scv = this->fvGeometry().scv(scvf.insideScvIdx());
246 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvf);
247
248 // set the Dirichlet conditions in residual and jacobian
249 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
250 {
251 static_assert(numEq == 1, "Not yet implemented for more than one vector-valued primary variable");
252 const int pvIdx = eqIdx;
253 const int componentIdx = scv.dofAxis();
254 if (bcTypes.isDirichlet(componentIdx))
255 applyDirichlet(scv, std::array<Scalar,1>{{dirichletValues[componentIdx]}}, eqIdx, pvIdx);
256 }
257 }
258 }
259 }
260 }
261 }
262
267 template<class... Args>
268 void maybeUpdateCouplingContext(Args&&...) {}
269
274 template<class... Args>
276};
277
287template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
289
295template<class TypeTag, class Assembler, class Implementation>
296class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
298 TypeTag, Assembler,
299 Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
300 /*implicit=*/true
301>
302{
306 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
308 using FVElementGeometry = typename GridGeometry::LocalView;
309 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
310 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
314
315 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
317
318public:
319
320 using LocalResidual = typename ParentType::LocalResidual;
321 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
322 using ParentType::ParentType;
323
330 template <class PartialReassembler = DefaultPartialReassembler>
331 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
332 const PartialReassembler* partialReassembler = nullptr)
333 {
334 // get some aliases for convenience
335 const auto& problem = this->asImp_().problem();
336 const auto& element = this->element();
337 const auto& fvGeometry = this->fvGeometry();
338 const auto& curSol = this->asImp_().curSol();
339 auto&& curElemVolVars = this->curElemVolVars();
340
341 // get the vector of the actual element residuals
342 const auto origResiduals = this->evalLocalResidual();
343
345 // //
346 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
347 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
348 // actual element. In the actual element we evaluate the derivative of the entire residual. //
349 // //
351
352 // one residual per element facet
353 const auto numElementResiduals = fvGeometry.numScv();
354
355 // create the vector storing the partial derivatives
356 ElementResidualVector partialDerivs(numElementResiduals);
357
358 const auto evalSource = [&](ElementResidualVector& residual, const SubControlVolume& scv)
359 {
360 this->localResidual().evalSource(residual, problem, element, fvGeometry, curElemVolVars, scv);
361 };
362
363 const auto evalStorage = [&](ElementResidualVector& residual, const SubControlVolume& scv)
364 {
365 this->localResidual().evalStorage(residual, problem, element, fvGeometry, this->prevElemVolVars(), curElemVolVars, scv);
366 };
367
368 const auto evalFlux = [&](ElementResidualVector& residual, const SubControlVolumeFace& scvf)
369 {
370 if (!scvf.processorBoundary())
371 this->localResidual().evalFlux(residual, problem, element, fvGeometry, curElemVolVars, this->elemBcTypes(), this->elemFluxVarsCache(), scvf);
372 };
373
374 const auto evalDerivative = [&] (const auto& scvI, const auto& scvJ)
375 {
376 // derivative w.r.t. own DOF
377 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
378 {
379 partialDerivs = 0.0;
380 const auto& otherElement = fvGeometry.gridGeometry().element(scvJ.elementIndex());
381 auto otherElemSol = elementSolution(otherElement, curSol, fvGeometry.gridGeometry()); // TODO allow selective creation of elemsol (for one scv)
382 auto& curOtherVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scvJ);
383 const VolumeVariables origOtherVolVars(curOtherVolVars);
384
385 auto evalResiduals = [&](Scalar priVar)
386 {
387 // update the volume variables and compute element residual
388 otherElemSol[scvJ.localDofIndex()][pvIdx] = priVar;
389 curOtherVolVars.update(otherElemSol, problem, otherElement, scvJ);
390 this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx);
391
392 ElementResidualVector residual(numElementResiduals);
393 residual = 0;
394
395 evalSource(residual, scvI);
396
397 if (!this->assembler().isStationaryProblem())
398 evalStorage(residual, scvI);
399
400 for (const auto& scvf : scvfs(fvGeometry, scvI))
401 evalFlux(residual, scvf);
402
403 return residual;
404 };
405
406 // derive the residuals numerically
407 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
408 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
409 NumericDifferentiation::partialDerivative(evalResiduals, otherElemSol[scvJ.localDofIndex()][pvIdx], partialDerivs, origResiduals,
410 eps_(otherElemSol[scvJ.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
411
412 const auto updateJacobian = [&]()
413 {
414 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
415 {
416 // A[i][col][eqIdx][pvIdx] is the rate of change of
417 // the residual of equation 'eqIdx' at dof 'i'
418 // depending on the primary variable 'pvIdx' at dof
419 // 'col'.
420 A[scvI.dofIndex()][scvJ.dofIndex()][eqIdx][pvIdx] += partialDerivs[scvI.localDofIndex()][eqIdx];
421 }
422 };
423
424 if (element.partitionType() == Dune::InteriorEntity)
425 updateJacobian();
426 else
427 {
428 const auto localIdxI = scvI.indexInElement();
429 const auto& facetI = element.template subEntity <1> (localIdxI);
430
431 if (facetI.partitionType() == Dune::BorderEntity)
432 updateJacobian();
433 }
434
435 // restore the original state of the scv's volume variables
436 curOtherVolVars = origOtherVolVars;
437
438 // restore the original element solution
439 otherElemSol[scvJ.localDofIndex()][pvIdx] = curSol[scvJ.dofIndex()][pvIdx];
440 this->asImp_().maybeUpdateCouplingContext(scvJ, otherElemSol, pvIdx);
441 // TODO additional dof dependencies
442 }
443 };
444
445 // calculation of the derivatives
446 for (const auto& scvI : scvs(fvGeometry))
447 {
448 // derivative w.r.t. own DOFs
449 evalDerivative(scvI, scvI);
450
451 // derivative w.r.t. other DOFs
452 const auto& otherScvIndices = fvGeometry.gridGeometry().connectivityMap()[scvI.index()];
453 for (const auto globalJ : otherScvIndices)
454 evalDerivative(scvI, fvGeometry.scv(globalJ));
455 }
456
457 // evaluate additional derivatives that might arise from the coupling
458 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
459
460 return origResiduals;
461 }
462};
463
464
470template<class TypeTag, class Assembler, class Implementation>
471class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation>
473 TypeTag, Assembler,
474 Detail::NonVoidOrDefault_t<Implementation, FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
475 /*implicit=*/false
476>
477{
478 using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>;
479 using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>;
481 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
485
486 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
487
488public:
489 using LocalResidual = typename ParentType::LocalResidual;
490 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
491 using ParentType::ParentType;
492
499 template <class PartialReassembler = DefaultPartialReassembler>
500 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
501 const PartialReassembler* partialReassembler = nullptr)
502 {
503 if (partialReassembler)
504 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
505
506 // get some aliases for convenience
507 const auto& problem = this->asImp_().problem();
508 const auto& element = this->element();
509 const auto& fvGeometry = this->fvGeometry();
510 const auto& curSol = this->asImp_().curSol();
511 auto&& curElemVolVars = this->curElemVolVars();
512
513 // get the vector of the actual element residuals
514 const auto origResiduals = this->evalLocalResidual();
515 const auto origStorageResiduals = this->evalLocalStorageResidual();
516
518 // //
519 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
520 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
521 // actual element. In the actual element we evaluate the derivative of the entire residual. //
522 // //
524
525 // create the element solution
526 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
527
528 // create the vector storing the partial derivatives
529 ElementResidualVector partialDerivs(fvGeometry.numScv());
530
531 // calculation of the derivatives
532 for (const auto& scv : scvs(fvGeometry))
533 {
534 // dof index and corresponding actual pri vars
535 const auto dofIdx = scv.dofIndex();
536 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
537 const VolumeVariables origVolVars(curVolVars);
538
539 // calculate derivatives w.r.t to the privars at the dof at hand
540 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
541 {
542 partialDerivs = 0.0;
543
544 auto evalStorage = [&](Scalar priVar)
545 {
546 // auto partialDerivsTmp = partialDerivs;
547 elemSol[scv.localDofIndex()][pvIdx] = priVar;
548 curVolVars.update(elemSol, problem, element, scv);
549 return this->evalLocalStorageResidual();
550 };
551
552 // derive the residuals numerically
553 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
554 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
555 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
556 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
557
558 // update the global stiffness matrix with the current partial derivatives
559 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
560 {
561 // A[i][col][eqIdx][pvIdx] is the rate of change of
562 // the residual of equation 'eqIdx' at dof 'i'
563 // depending on the primary variable 'pvIdx' at dof
564 // 'col'.
565 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
566 }
567
568 // restore the original state of the scv's volume variables
569 curVolVars = origVolVars;
570
571 // restore the original element solution
572 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
573 }
574 }
575 return origResiduals;
576 }
577};
578
584template<class TypeTag, class Assembler, class Implementation>
585class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation>
587 TypeTag, Assembler,
588 FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>,
589 /*implicit=*/true
590>
591{
592 using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>;
593 using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, true>;
597
598 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
599
600public:
601 using LocalResidual = typename ParentType::LocalResidual;
602 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
603 using ParentType::ParentType;
604
611 template <class PartialReassembler = DefaultPartialReassembler>
612 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
613 const PartialReassembler* partialReassembler = nullptr)
614 {
615 if (partialReassembler)
616 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
617
618 // get some aliases for convenience
619 const auto& element = this->element();
620 const auto& fvGeometry = this->fvGeometry();
621 const auto& problem = this->asImp_().problem();
622 const auto& curElemVolVars = this->curElemVolVars();
623 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
624
625 // get the vecor of the actual element residuals
626 const auto origResiduals = this->evalLocalResidual();
627
629 // //
630 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
631 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
632 // actual element. In the actual element we evaluate the derivative of the entire residual. //
633 // //
635
636 // calculation of the source and storage derivatives
637 for (const auto& scv : scvs(fvGeometry))
638 {
639 // dof index and corresponding actual pri vars
640 const auto dofIdx = scv.dofIndex();
641 const auto& volVars = curElemVolVars[scv];
642
643 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
644 // only if the problem is instationary we add derivative of storage term
645 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
646 if (!this->assembler().isStationaryProblem())
647 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
648 problem,
649 element,
650 fvGeometry,
651 volVars,
652 scv);
653
654 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
655 // add source term derivatives
656 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
657 problem,
658 element,
659 fvGeometry,
660 volVars,
661 scv);
662 }
663
664 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
665 for (const auto& scvf : scvfs(fvGeometry))
666 {
667 if (!scvf.boundary())
668 {
669 // add flux term derivatives
670 this->localResidual().addFluxDerivatives(A,
671 problem,
672 element,
673 fvGeometry,
674 curElemVolVars,
675 elemFluxVarsCache,
676 scvf);
677 }
678
679 // the boundary gets special treatment to simplify
680 // for the user
681 else
682 {
683 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
684 if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
685 {
686 // add flux term derivatives
687 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
688 problem,
689 element,
690 fvGeometry,
691 curElemVolVars,
692 elemFluxVarsCache,
693 scvf);
694 }
695 }
696 }
697
698 return origResiduals;
699 }
700};
701
707template<class TypeTag, class Assembler, class Implementation>
708class FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation>
710 TypeTag, Assembler,
711 FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>,
712 /*implicit=*/false
713>
714{
715 using ThisType = FaceCenteredLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>;
716 using ParentType = FaceCenteredLocalAssemblerBase<TypeTag, Assembler, Detail::NonVoidOrDefault_t<Implementation, ThisType>, false>;
720
721 enum { numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq() };
722
723public:
724 using LocalResidual = typename ParentType::LocalResidual;
725 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
726 using ParentType::ParentType;
727
734 template <class PartialReassembler = DefaultPartialReassembler>
735 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
736 const PartialReassembler* partialReassembler = nullptr)
737 {
738 if (partialReassembler)
739 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
740
741 // get some aliases for convenience
742 const auto& element = this->element();
743 const auto& fvGeometry = this->fvGeometry();
744 const auto& problem = this->asImp_().problem();
745 const auto& curElemVolVars = this->curElemVolVars();
746
747 // get the vector of the actual element residuals
748 const auto origResiduals = this->evalLocalResidual();
749
751 // //
752 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
753 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
754 // actual element. In the actual element we evaluate the derivative of the entire residual. //
755 // //
757
758 // calculation of the source and storage derivatives
759 for (const auto& scv : scvs(fvGeometry))
760 {
761 // dof index and corresponding actual pri vars
762 const auto dofIdx = scv.dofIndex();
763 const auto& volVars = curElemVolVars[scv];
764
765 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
766 // only if the problem is instationary we add derivative of storage term
767 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
768 problem,
769 element,
770 fvGeometry,
771 volVars,
772 scv);
773 }
774
775 return origResiduals;
776 }
777};
778
779} // end namespace Dumux
780
781#endif
A base class for all local assemblers.
void bindLocalViews()
Definition assembly/fvlocalassemblerbase.hh:173
ElementVolumeVariables & curElemVolVars()
Definition assembly/fvlocalassemblerbase.hh:253
ElementBoundaryTypes & elemBcTypes()
Definition assembly/fvlocalassemblerbase.hh:269
Implementation & asImp_()
Definition assembly/fvlocalassemblerbase.hh:297
ElementResidualVector evalLocalResidual() const
Definition assembly/fvlocalassemblerbase.hh:108
const Problem & problem() const
Definition assembly/fvlocalassemblerbase.hh:229
FVLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol)
Definition assembly/fvlocalassemblerbase.hh:61
FVElementGeometry & fvGeometry()
Definition assembly/fvlocalassemblerbase.hh:249
bool elementIsGhost() const
Definition assembly/fvlocalassemblerbase.hh:241
const Element & element() const
Definition assembly/fvlocalassemblerbase.hh:237
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:79
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition fclocalassembler.hh:220
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:173
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition fclocalassembler.hh:268
void bindLocalViews()
Definition fclocalassembler.hh:68
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition fclocalassembler.hh:232
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition fclocalassembler.hh:275
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition fclocalassembler.hh:197
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition fclocalassembler.hh:288
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:50
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.
The global face variables class for staggered models.
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition diffmethod.hh:25
@ analytic
Definition diffmethod.hh:26
@ numeric
Definition diffmethod.hh:26
@ green
does not need to be reassembled
Definition entitycolor.hh:40
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
@ element
Definition fieldtype.hh:23
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
Distance implementation details.
Definition cvfelocalresidual.hh:25
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