version 3.10-dev
assembly/cvfelocalassembler.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_CVFE_LOCAL_ASSEMBLER_HH
14#define DUMUX_CVFE_LOCAL_ASSEMBLER_HH
15
16#include <dune/common/exceptions.hh>
17#include <dune/common/hybridutilities.hh>
18#include <dune/common/reservedvector.hh>
19#include <dune/grid/common/gridenums.hh>
20#include <dune/istl/matrixindexset.hh>
21#include <dune/istl/bvector.hh>
22
26
32
34
35#include "volvardeflectionhelper_.hh"
36
37namespace Dumux {
38
39#ifndef DOXYGEN
40namespace Detail::CVFE {
41
42struct NoOperator
43{
44 template<class... Args>
45 constexpr void operator()(Args&&...) const {}
46};
47
48template<class X, class Y>
49using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>;
50
51} // end namespace Detail
52#endif // DOXYGEN
53
63template<class TypeTag, class Assembler, class Implementation, bool implicit>
64class CVFELocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
65{
69 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
70 using SolutionVector = typename Assembler::SolutionVector;
71
72 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
74
75public:
76
77 using ParentType::ParentType;
78
80 {
82 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
83 }
84
85
90 template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::CVFE::NoOperator>
91 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
92 const PartialReassembler* partialReassembler = nullptr,
93 const CouplingFunction& maybeAssembleCouplingBlocks = {})
94 {
95 this->asImp_().bindLocalViews();
96 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
97 if (partialReassembler
98 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
99 {
100 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
101 for (const auto& scv : scvs(this->fvGeometry()))
102 res[scv.dofIndex()] += residual[scv.localDofIndex()];
103
104 // assemble the coupling blocks for coupled models (does nothing if not coupled)
105 maybeAssembleCouplingBlocks(residual);
106 }
107 else if (!this->elementIsGhost())
108 {
109 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
110 for (const auto& scv : scvs(this->fvGeometry()))
111 res[scv.dofIndex()] += residual[scv.localDofIndex()];
112
113 // assemble the coupling blocks for coupled models (does nothing if not coupled)
114 maybeAssembleCouplingBlocks(residual);
115 }
116 else
117 {
118 // Treatment of ghost elements
119 assert(this->elementIsGhost());
120
121 // handle dofs per codimension
122 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
123 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim>{}, [&](auto d)
124 {
125 constexpr int codim = dim - d;
126 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
127 for (int idx = 0; idx < localCoeffs.size(); ++idx)
128 {
129 const auto& localKey = localCoeffs.localKey(idx);
130
131 // skip if we are not handling this codim right now
132 if (localKey.codim() != codim)
133 continue;
134
135 // do not change the non-ghost entities
136 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
137 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
138 continue;
139
140 // WARNING: this only works if the mapping from codim+subEntity to
141 // global dofIndex is unique (on dof per entity of this codim).
142 // For more general mappings, we should use a proper local-global mapping here.
143 // For example through dune-functions.
144 const auto dofIndex = gridGeometry.dofMapper().index(entity);
145
146 // this might be a vector-valued dof
147 using BlockType = typename JacobianMatrix::block_type;
148 BlockType &J = jac[dofIndex][dofIndex];
149 for (int j = 0; j < BlockType::rows; ++j)
150 J[j][j] = 1.0;
151
152 // set residual for the ghost dof
153 res[dofIndex] = 0;
154 }
155 });
156 }
157
158 auto applyDirichlet = [&] (const auto& scvI,
159 const auto& dirichletValues,
160 const auto eqIdx,
161 const auto pvIdx)
162 {
163 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
164
165 auto& row = jac[scvI.dofIndex()];
166 for (auto col = row.begin(); col != row.end(); ++col)
167 row[col.index()][eqIdx] = 0.0;
168
169 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
170
171 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
172 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
173 {
174 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
175 res[periodicDof][eqIdx] = this->asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx];
176
177 auto& rowP = jac[periodicDof];
178 for (auto col = rowP.begin(); col != rowP.end(); ++col)
179 rowP[col.index()][eqIdx] = 0.0;
180
181 rowP[periodicDof][eqIdx][pvIdx] = 1.0;
182 }
183 };
184
185 this->asImp_().enforceDirichletConstraints(applyDirichlet);
186 }
187
192 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
193 {
194 this->asImp_().bindLocalViews();
195 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
196
197 auto applyDirichlet = [&] (const auto& scvI,
198 const auto& dirichletValues,
199 const auto eqIdx,
200 const auto pvIdx)
201 {
202 auto& row = jac[scvI.dofIndex()];
203 for (auto col = row.begin(); col != row.end(); ++col)
204 row[col.index()][eqIdx] = 0.0;
205
206 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
207 };
208
209 this->asImp_().enforceDirichletConstraints(applyDirichlet);
210 }
211
215 template <class ResidualVector>
216 void assembleResidual(ResidualVector& res)
217 {
218 this->asImp_().bindLocalViews();
219 const auto residual = this->evalLocalResidual();
220
221 for (const auto& scv : scvs(this->fvGeometry()))
222 res[scv.dofIndex()] += residual[scv.localDofIndex()];
223
224 auto applyDirichlet = [&] (const auto& scvI,
225 const auto& dirichletValues,
226 const auto eqIdx,
227 const auto pvIdx)
228 {
229 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
230 };
231
232 this->asImp_().enforceDirichletConstraints(applyDirichlet);
233 }
234
236 template<typename ApplyFunction>
237 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
238 {
239 // enforce Dirichlet boundary conditions
240 this->asImp_().evalDirichletBoundaries(applyDirichlet);
241 // take care of internal Dirichlet constraints (if enabled)
242 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
243 }
244
248 template< typename ApplyDirichletFunctionType >
249 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
250 {
251 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
252 // and set the residual to (privar - dirichletvalue)
253 if (this->elemBcTypes().hasDirichlet())
254 {
255 for (const auto& scvI : scvs(this->fvGeometry()))
256 {
257 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
258 if (bcTypes.hasDirichlet())
259 {
260 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
261
262 // set the Dirichlet conditions in residual and jacobian
263 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
264 {
265 if (bcTypes.isDirichlet(eqIdx))
266 {
267 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
268 assert(0 <= pvIdx && pvIdx < numEq);
269 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
270 }
271 }
272 }
273 }
274 }
275 }
276
281 template<class... Args>
282 void maybeUpdateCouplingContext(Args&&...) {}
283
288 template<class... Args>
290};
291
301template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
303
309template<class TypeTag, class Assembler, class Implementation>
310class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
311: public CVFELocalAssemblerBase<TypeTag, Assembler,
312 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
313 true>
314{
321
322 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
324
325 static constexpr bool enableGridFluxVarsCache
326 = GridVariables::GridFluxVariablesCache::cachingEnabled;
327 static constexpr bool solutionDependentFluxVarsCache
328 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
329
330public:
331
333 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
335
342 template <class PartialReassembler = DefaultPartialReassembler>
343 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
344 const PartialReassembler* partialReassembler = nullptr)
345 {
346 // get some aliases for convenience
347 const auto& element = this->element();
348 const auto& fvGeometry = this->fvGeometry();
349 const auto& curSol = this->asImp_().curSol();
350
351 auto&& curElemVolVars = this->curElemVolVars();
352 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
353
354 // get the vector of the actual element residuals
355 const auto origResiduals = this->evalLocalResidual();
356
358 // //
359 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
360 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
361 // actual element. In the actual element we evaluate the derivative of the entire residual. //
362 // //
364
365 // if all volvars in the stencil have to be updated or if it's enough to only update the
366 // volVars for the scv whose associated dof has been deflected
367 static const bool updateAllVolVars = getParamFromGroup<bool>(
368 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
369 );
370
371 // create the element solution
372 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
373
374 // create the vector storing the partial derivatives
375 ElementResidualVector partialDerivs(fvGeometry.numScv());
376
377 Detail::VolVarsDeflectionHelper deflectionHelper(
378 [&] (const auto& scv) -> VolumeVariables& {
379 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
380 },
381 fvGeometry,
382 updateAllVolVars
383 );
384
385 // calculation of the derivatives
386 for (const auto& scv : scvs(fvGeometry))
387 {
388 // dof index and corresponding actual pri vars
389 const auto dofIdx = scv.dofIndex();
390 deflectionHelper.setCurrent(scv);
391
392 // calculate derivatives w.r.t to the privars at the dof at hand
393 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
394 {
395 partialDerivs = 0.0;
396
397 auto evalResiduals = [&](Scalar priVar)
398 {
399 // update the volume variables and compute element residual
400 elemSol[scv.localDofIndex()][pvIdx] = priVar;
401 deflectionHelper.deflect(elemSol, scv, this->asImp_().problem());
402 if constexpr (solutionDependentFluxVarsCache)
403 {
404 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
405 if constexpr (enableGridFluxVarsCache)
406 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
407 }
408 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
409 return this->evalLocalResidual();
410 };
411
412 // derive the residuals numerically
413 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
414 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
415 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
416 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
417
418 // update the global stiffness matrix with the current partial derivatives
419 for (const auto& scvJ : scvs(fvGeometry))
420 {
421 // don't add derivatives for green dofs
422 if (!partialReassembler
423 || partialReassembler->dofColor(scvJ.dofIndex()) != EntityColor::green)
424 {
425 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
426 {
427 // A[i][col][eqIdx][pvIdx] is the rate of change of
428 // the residual of equation 'eqIdx' at dof 'i'
429 // depending on the primary variable 'pvIdx' at dof
430 // 'col'.
431 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
432 }
433 }
434 }
435
436 // restore the original state of the scv's volume variables
437 deflectionHelper.restore(scv);
438
439 // restore the original element solution
440 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
441 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
442 }
443 }
444
445 // restore original state of the flux vars cache in case of global caching.
446 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
447 if constexpr (enableGridFluxVarsCache)
448 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
449
450 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
451 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
452
453 return origResiduals;
454 }
455
456}; // implicit CVFEAssembler with numeric Jacobian
457
463template<class TypeTag, class Assembler, class Implementation>
464class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation>
465: public CVFELocalAssemblerBase<TypeTag, Assembler,
466 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
467 false>
468{
475
476 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
478
479public:
480
482 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
484
491 template <class PartialReassembler = DefaultPartialReassembler>
492 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
493 const PartialReassembler* partialReassembler = nullptr)
494 {
495 if (partialReassembler)
496 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
497
498 // get some aliases for convenience
499 const auto& element = this->element();
500 const auto& fvGeometry = this->fvGeometry();
501 const auto& curSol = this->asImp_().curSol();
502 auto&& curElemVolVars = this->curElemVolVars();
503
504 // get the vecor of the actual element residuals
505 const auto origResiduals = this->evalLocalResidual();
506 const auto origStorageResiduals = this->evalLocalStorageResidual();
507
509 // //
510 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
511 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
512 // actual element. In the actual element we evaluate the derivative of the entire residual. //
513 // //
515
516 // create the element solution
517 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
518
519 // create the vector storing the partial derivatives
520 ElementResidualVector partialDerivs(fvGeometry.numScv());
521
522 // calculation of the derivatives
523 for (const auto& scv : scvs(fvGeometry))
524 {
525 // dof index and corresponding actual pri vars
526 const auto dofIdx = scv.dofIndex();
527 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
528 const VolumeVariables origVolVars(curVolVars);
529
530 // calculate derivatives w.r.t to the privars at the dof at hand
531 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
532 {
533 partialDerivs = 0.0;
534
535 auto evalStorage = [&](Scalar priVar)
536 {
537 // auto partialDerivsTmp = partialDerivs;
538 elemSol[scv.localDofIndex()][pvIdx] = priVar;
539 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
540 return this->evalLocalStorageResidual();
541 };
542
543 // derive the residuals numerically
544 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
545 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
546 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
547 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
548
549 // update the global stiffness matrix with the current partial derivatives
550 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
551 {
552 // A[i][col][eqIdx][pvIdx] is the rate of change of
553 // the residual of equation 'eqIdx' at dof 'i'
554 // depending on the primary variable 'pvIdx' at dof
555 // 'col'.
556 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
557 }
558
559 // restore the original state of the scv's volume variables
560 curVolVars = origVolVars;
561
562 // restore the original element solution
563 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
564 // TODO additional dof dependencies
565 }
566 }
567
568 return origResiduals;
569 }
570}; // explicit CVFEAssembler with numeric Jacobian
571
577template<class TypeTag, class Assembler, class Implementation>
578class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation>
579: public CVFELocalAssemblerBase<TypeTag, Assembler,
580 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
581 true>
582{
587
588public:
589
591 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
593
600 template <class PartialReassembler = DefaultPartialReassembler>
601 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
602 const PartialReassembler* partialReassembler = nullptr)
603 {
604 if (partialReassembler)
605 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
606
607 // get some aliases for convenience
608 const auto& element = this->element();
609 const auto& fvGeometry = this->fvGeometry();
610 const auto& problem = this->asImp_().problem();
611 const auto& curElemVolVars = this->curElemVolVars();
612 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
613
614 // get the vecor of the actual element residuals
615 const auto origResiduals = this->evalLocalResidual();
616
618 // //
619 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
620 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
621 // actual element. In the actual element we evaluate the derivative of the entire residual. //
622 // //
624
625 // calculation of the source and storage derivatives
626 for (const auto& scv : scvs(fvGeometry))
627 {
628 // dof index and corresponding actual pri vars
629 const auto dofIdx = scv.dofIndex();
630 const auto& volVars = curElemVolVars[scv];
631
632 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
633 // only if the problem is instationary we add derivative of storage term
634 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
635 if (!this->assembler().isStationaryProblem())
636 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
637 problem,
638 element,
639 fvGeometry,
640 volVars,
641 scv);
642
643 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
644 // add source term derivatives
645 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
646 problem,
647 element,
648 fvGeometry,
649 volVars,
650 scv);
651 }
652
653 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
654 for (const auto& scvf : scvfs(fvGeometry))
655 {
656 if (!scvf.boundary())
657 {
658 // add flux term derivatives
659 this->localResidual().addFluxDerivatives(A,
660 problem,
661 element,
662 fvGeometry,
663 curElemVolVars,
664 elemFluxVarsCache,
665 scvf);
666 }
667
668 // the boundary gets special treatment to simplify
669 // for the user
670 else
671 {
672 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
673 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
674 {
675 // add flux term derivatives
676 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
677 problem,
678 element,
679 fvGeometry,
680 curElemVolVars,
681 elemFluxVarsCache,
682 scvf);
683 }
684 }
685 }
686
687 return origResiduals;
688 }
689
690}; // implicit CVFEAssembler with analytic Jacobian
691
697template<class TypeTag, class Assembler, class Implementation>
698class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation>
699: public CVFELocalAssemblerBase<TypeTag, Assembler,
700 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
701 false>
702{
707
708public:
709
711 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
713
720 template <class PartialReassembler = DefaultPartialReassembler>
721 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
722 const PartialReassembler* partialReassembler = nullptr)
723 {
724 if (partialReassembler)
725 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
726
727 // get some aliases for convenience
728 const auto& element = this->element();
729 const auto& fvGeometry = this->fvGeometry();
730 const auto& problem = this->asImp_().problem();
731 const auto& curElemVolVars = this->curElemVolVars();
732
733 // get the vecor of the actual element residuals
734 const auto origResiduals = this->evalLocalResidual();
735
737 // //
738 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
739 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
740 // actual element. In the actual element we evaluate the derivative of the entire residual. //
741 // //
743
744 // calculation of the source and storage derivatives
745 for (const auto& scv : scvs(fvGeometry))
746 {
747 // dof index and corresponding actual pri vars
748 const auto dofIdx = scv.dofIndex();
749 const auto& volVars = curElemVolVars[scv];
750
751 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
752 // only if the problem is instationary we add derivative of storage term
753 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
754 problem,
755 element,
756 fvGeometry,
757 volVars,
758 scv);
759 }
760
761 return origResiduals;
762 }
763
764}; // explicit CVFEAssembler with analytic Jacobian
765
766} // end namespace Dumux
767
768#endif
Control volume finite element local assembler using analytic differentiation and explicit time discre...
Definition: assembly/cvfelocalassembler.hh:702
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:711
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: assembly/cvfelocalassembler.hh:721
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:710
Control volume finite element local assembler using analytic differentiation and implicit time discre...
Definition: assembly/cvfelocalassembler.hh:582
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: assembly/cvfelocalassembler.hh:601
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:591
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:590
Control volume finite element local assembler using numeric differentiation and implicit time discret...
Definition: assembly/cvfelocalassembler.hh:314
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:333
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:332
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: assembly/cvfelocalassembler.hh:343
Control volume finite element local assembler using numeric differentiation and explicit time discret...
Definition: assembly/cvfelocalassembler.hh:468
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:482
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: assembly/cvfelocalassembler.hh:492
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:481
A base class for all local CVFE assemblers.
Definition: assembly/cvfelocalassembler.hh:65
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: assembly/cvfelocalassembler.hh:237
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: assembly/cvfelocalassembler.hh:216
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: assembly/cvfelocalassembler.hh:282
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cvfelocalassembler.hh:192
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: assembly/cvfelocalassembler.hh:249
void assembleJacobianAndResidual(JacobianMatrix &jac, ResidualVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr, const CouplingFunction &maybeAssembleCouplingBlocks={})
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: assembly/cvfelocalassembler.hh:91
void bindLocalViews()
Definition: assembly/cvfelocalassembler.hh:79
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: assembly/cvfelocalassembler.hh:289
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:302
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
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
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
Defines all properties used in Dumux.
The local element solution class for control-volume finite element methods.
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
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.