version 3.8
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->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
176 const auto end = jac[periodicDof].end();
177 for (auto it = jac[periodicDof].begin(); it != end; ++it)
178 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
179 }
180 };
181
182 this->asImp_().enforceDirichletConstraints(applyDirichlet);
183 }
184
189 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
190 {
191 this->asImp_().bindLocalViews();
192 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
193
194 auto applyDirichlet = [&] (const auto& scvI,
195 const auto& dirichletValues,
196 const auto eqIdx,
197 const auto pvIdx)
198 {
199 auto& row = jac[scvI.dofIndex()];
200 for (auto col = row.begin(); col != row.end(); ++col)
201 row[col.index()][eqIdx] = 0.0;
202
203 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
204 };
205
206 this->asImp_().enforceDirichletConstraints(applyDirichlet);
207 }
208
212 template <class ResidualVector>
213 void assembleResidual(ResidualVector& res)
214 {
215 this->asImp_().bindLocalViews();
216 const auto residual = this->evalLocalResidual();
217
218 for (const auto& scv : scvs(this->fvGeometry()))
219 res[scv.dofIndex()] += residual[scv.localDofIndex()];
220
221 auto applyDirichlet = [&] (const auto& scvI,
222 const auto& dirichletValues,
223 const auto eqIdx,
224 const auto pvIdx)
225 {
226 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
227 };
228
229 this->asImp_().enforceDirichletConstraints(applyDirichlet);
230 }
231
233 template<typename ApplyFunction>
234 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
235 {
236 // enforce Dirichlet boundary conditions
237 this->asImp_().evalDirichletBoundaries(applyDirichlet);
238 // take care of internal Dirichlet constraints (if enabled)
239 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
240 }
241
245 template< typename ApplyDirichletFunctionType >
246 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
247 {
248 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
249 // and set the residual to (privar - dirichletvalue)
250 if (this->elemBcTypes().hasDirichlet())
251 {
252 for (const auto& scvI : scvs(this->fvGeometry()))
253 {
254 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
255 if (bcTypes.hasDirichlet())
256 {
257 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
258
259 // set the Dirichlet conditions in residual and jacobian
260 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
261 {
262 if (bcTypes.isDirichlet(eqIdx))
263 {
264 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
265 assert(0 <= pvIdx && pvIdx < numEq);
266 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
267 }
268 }
269 }
270 }
271 }
272 }
273
278 template<class... Args>
279 void maybeUpdateCouplingContext(Args&&...) {}
280
285 template<class... Args>
287};
288
298template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
300
306template<class TypeTag, class Assembler, class Implementation>
307class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
308: public CVFELocalAssemblerBase<TypeTag, Assembler,
309 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
310 true>
311{
318
319 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
321
322 static constexpr bool enableGridFluxVarsCache
323 = GridVariables::GridFluxVariablesCache::cachingEnabled;
324 static constexpr bool solutionDependentFluxVarsCache
325 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
326
327public:
328
330 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
332
339 template <class PartialReassembler = DefaultPartialReassembler>
340 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
341 const PartialReassembler* partialReassembler = nullptr)
342 {
343 // get some aliases for convenience
344 const auto& element = this->element();
345 const auto& fvGeometry = this->fvGeometry();
346 const auto& curSol = this->asImp_().curSol();
347
348 auto&& curElemVolVars = this->curElemVolVars();
349 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
350
351 // get the vector of the actual element residuals
352 const auto origResiduals = this->evalLocalResidual();
353
355 // //
356 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
357 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
358 // actual element. In the actual element we evaluate the derivative of the entire residual. //
359 // //
361
362 // if all volvars in the stencil have to be updated or if it's enough to only update the
363 // volVars for the scv whose associated dof has been deflected
364 static const bool updateAllVolVars = getParamFromGroup<bool>(
365 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
366 );
367
368 // create the element solution
369 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
370
371 // create the vector storing the partial derivatives
372 ElementResidualVector partialDerivs(fvGeometry.numScv());
373
374 Detail::VolVarsDeflectionHelper deflectionHelper(
375 [&] (const auto& scv) -> VolumeVariables& {
376 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
377 },
378 fvGeometry,
379 updateAllVolVars
380 );
381
382 // calculation of the derivatives
383 for (const auto& scv : scvs(fvGeometry))
384 {
385 // dof index and corresponding actual pri vars
386 const auto dofIdx = scv.dofIndex();
387 deflectionHelper.setCurrent(scv);
388
389 // calculate derivatives w.r.t to the privars at the dof at hand
390 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
391 {
392 partialDerivs = 0.0;
393
394 auto evalResiduals = [&](Scalar priVar)
395 {
396 // update the volume variables and compute element residual
397 elemSol[scv.localDofIndex()][pvIdx] = priVar;
398 deflectionHelper.deflect(elemSol, scv, this->asImp_().problem());
399 if constexpr (solutionDependentFluxVarsCache)
400 {
401 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
402 if constexpr (enableGridFluxVarsCache)
403 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
404 }
405 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
406 return this->evalLocalResidual();
407 };
408
409 // derive the residuals numerically
410 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
411 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
412 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
413 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
414
415 // update the global stiffness matrix with the current partial derivatives
416 for (const auto& scvJ : scvs(fvGeometry))
417 {
418 // don't add derivatives for green dofs
419 if (!partialReassembler
420 || partialReassembler->dofColor(scvJ.dofIndex()) != EntityColor::green)
421 {
422 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
423 {
424 // A[i][col][eqIdx][pvIdx] is the rate of change of
425 // the residual of equation 'eqIdx' at dof 'i'
426 // depending on the primary variable 'pvIdx' at dof
427 // 'col'.
428 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
429 }
430 }
431 }
432
433 // restore the original state of the scv's volume variables
434 deflectionHelper.restore(scv);
435
436 // restore the original element solution
437 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
438 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
439 }
440 }
441
442 // restore original state of the flux vars cache in case of global caching.
443 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
444 if constexpr (enableGridFluxVarsCache)
445 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
446
447 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
448 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
449
450 return origResiduals;
451 }
452
453}; // implicit CVFEAssembler with numeric Jacobian
454
460template<class TypeTag, class Assembler, class Implementation>
461class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation>
462: public CVFELocalAssemblerBase<TypeTag, Assembler,
463 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
464 false>
465{
472
473 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
475
476public:
477
479 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
481
488 template <class PartialReassembler = DefaultPartialReassembler>
489 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
490 const PartialReassembler* partialReassembler = nullptr)
491 {
492 if (partialReassembler)
493 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
494
495 // get some aliases for convenience
496 const auto& element = this->element();
497 const auto& fvGeometry = this->fvGeometry();
498 const auto& curSol = this->asImp_().curSol();
499 auto&& curElemVolVars = this->curElemVolVars();
500
501 // get the vecor of the actual element residuals
502 const auto origResiduals = this->evalLocalResidual();
503 const auto origStorageResiduals = this->evalLocalStorageResidual();
504
506 // //
507 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
508 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
509 // actual element. In the actual element we evaluate the derivative of the entire residual. //
510 // //
512
513 // create the element solution
514 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
515
516 // create the vector storing the partial derivatives
517 ElementResidualVector partialDerivs(fvGeometry.numScv());
518
519 // calculation of the derivatives
520 for (const auto& scv : scvs(fvGeometry))
521 {
522 // dof index and corresponding actual pri vars
523 const auto dofIdx = scv.dofIndex();
524 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
525 const VolumeVariables origVolVars(curVolVars);
526
527 // calculate derivatives w.r.t to the privars at the dof at hand
528 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
529 {
530 partialDerivs = 0.0;
531
532 auto evalStorage = [&](Scalar priVar)
533 {
534 // auto partialDerivsTmp = partialDerivs;
535 elemSol[scv.localDofIndex()][pvIdx] = priVar;
536 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
537 return this->evalLocalStorageResidual();
538 };
539
540 // derive the residuals numerically
541 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
542 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
543 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
544 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
545
546 // update the global stiffness matrix with the current partial derivatives
547 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
548 {
549 // A[i][col][eqIdx][pvIdx] is the rate of change of
550 // the residual of equation 'eqIdx' at dof 'i'
551 // depending on the primary variable 'pvIdx' at dof
552 // 'col'.
553 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
554 }
555
556 // restore the original state of the scv's volume variables
557 curVolVars = origVolVars;
558
559 // restore the original element solution
560 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
561 // TODO additional dof dependencies
562 }
563 }
564
565 return origResiduals;
566 }
567}; // explicit CVFEAssembler with numeric Jacobian
568
574template<class TypeTag, class Assembler, class Implementation>
575class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation>
576: public CVFELocalAssemblerBase<TypeTag, Assembler,
577 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
578 true>
579{
584
585public:
586
588 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
590
597 template <class PartialReassembler = DefaultPartialReassembler>
598 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
599 const PartialReassembler* partialReassembler = nullptr)
600 {
601 if (partialReassembler)
602 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
603
604 // get some aliases for convenience
605 const auto& element = this->element();
606 const auto& fvGeometry = this->fvGeometry();
607 const auto& problem = this->asImp_().problem();
608 const auto& curElemVolVars = this->curElemVolVars();
609 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
610
611 // get the vecor of the actual element residuals
612 const auto origResiduals = this->evalLocalResidual();
613
615 // //
616 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
617 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
618 // actual element. In the actual element we evaluate the derivative of the entire residual. //
619 // //
621
622 // calculation of the source and storage derivatives
623 for (const auto& scv : scvs(fvGeometry))
624 {
625 // dof index and corresponding actual pri vars
626 const auto dofIdx = scv.dofIndex();
627 const auto& volVars = curElemVolVars[scv];
628
629 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
630 // only if the problem is instationary we add derivative of storage term
631 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
632 if (!this->assembler().isStationaryProblem())
633 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
634 problem,
635 element,
636 fvGeometry,
637 volVars,
638 scv);
639
640 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
641 // add source term derivatives
642 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
643 problem,
644 element,
645 fvGeometry,
646 volVars,
647 scv);
648 }
649
650 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
651 for (const auto& scvf : scvfs(fvGeometry))
652 {
653 if (!scvf.boundary())
654 {
655 // add flux term derivatives
656 this->localResidual().addFluxDerivatives(A,
657 problem,
658 element,
659 fvGeometry,
660 curElemVolVars,
661 elemFluxVarsCache,
662 scvf);
663 }
664
665 // the boundary gets special treatment to simplify
666 // for the user
667 else
668 {
669 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
670 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
671 {
672 // add flux term derivatives
673 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
674 problem,
675 element,
676 fvGeometry,
677 curElemVolVars,
678 elemFluxVarsCache,
679 scvf);
680 }
681 }
682 }
683
684 return origResiduals;
685 }
686
687}; // implicit CVFEAssembler with analytic Jacobian
688
694template<class TypeTag, class Assembler, class Implementation>
695class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation>
696: public CVFELocalAssemblerBase<TypeTag, Assembler,
697 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
698 false>
699{
704
705public:
706
708 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
710
717 template <class PartialReassembler = DefaultPartialReassembler>
718 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
719 const PartialReassembler* partialReassembler = nullptr)
720 {
721 if (partialReassembler)
722 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
723
724 // get some aliases for convenience
725 const auto& element = this->element();
726 const auto& fvGeometry = this->fvGeometry();
727 const auto& problem = this->asImp_().problem();
728 const auto& curElemVolVars = this->curElemVolVars();
729
730 // get the vecor of the actual element residuals
731 const auto origResiduals = this->evalLocalResidual();
732
734 // //
735 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
736 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
737 // actual element. In the actual element we evaluate the derivative of the entire residual. //
738 // //
740
741 // calculation of the source and storage derivatives
742 for (const auto& scv : scvs(fvGeometry))
743 {
744 // dof index and corresponding actual pri vars
745 const auto dofIdx = scv.dofIndex();
746 const auto& volVars = curElemVolVars[scv];
747
748 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
749 // only if the problem is instationary we add derivative of storage term
750 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
751 problem,
752 element,
753 fvGeometry,
754 volVars,
755 scv);
756 }
757
758 return origResiduals;
759 }
760
761}; // explicit CVFEAssembler with analytic Jacobian
762
763} // end namespace Dumux
764
765#endif
Control volume finite element local assembler using analytic differentiation and explicit time discre...
Definition: assembly/cvfelocalassembler.hh:699
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:708
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:718
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:707
Control volume finite element local assembler using analytic differentiation and implicit time discre...
Definition: assembly/cvfelocalassembler.hh:579
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:598
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:588
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:587
Control volume finite element local assembler using numeric differentiation and implicit time discret...
Definition: assembly/cvfelocalassembler.hh:311
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:330
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:329
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cvfelocalassembler.hh:340
Control volume finite element local assembler using numeric differentiation and explicit time discret...
Definition: assembly/cvfelocalassembler.hh:465
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: assembly/cvfelocalassembler.hh:479
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:489
typename ParentType::LocalResidual LocalResidual
Definition: assembly/cvfelocalassembler.hh:478
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:234
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: assembly/cvfelocalassembler.hh:213
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: assembly/cvfelocalassembler.hh:279
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:189
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: assembly/cvfelocalassembler.hh:246
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:286
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:299
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: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
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.