version 3.11-dev
Loading...
Searching...
No Matches
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-FileCopyrightText: 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
27#include <dumux/common/typetraits/localdofs_.hh>
28#include <dumux/common/typetraits/boundary_.hh>
29
35
38
39#include "cvfevolvarsdeflectionpolicy_.hh"
40
41namespace Dumux {
42
43#ifndef DOXYGEN
44namespace Detail::CVFE {
45
46struct NoOperator
47{
48 template<class... Args>
49 constexpr void operator()(Args&&...) const {}
50};
51
52template<class X, class Y>
53using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>;
54
55} // end namespace Detail
56#endif // DOXYGEN
57
67template<class TypeTag, class Assembler, class Implementation, bool implicit>
68class CVFELocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
69{
74 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
75 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
77 using SolutionVector = typename Assembler::SolutionVector;
79 using FVElementGeometry = typename GridGeometry::LocalView;
80
81 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
82 static constexpr int dim = GridGeometry::GridView::dimension;
83
84public:
85
86 using ParentType::ParentType;
87
89 {
91 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
92 }
93
94
99 template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::CVFE::NoOperator>
100 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
101 const PartialReassembler* partialReassembler = nullptr,
102 const CouplingFunction& maybeAssembleCouplingBlocks = {})
103 {
104 this->asImp_().bindLocalViews();
105 const auto eIdxGlobal = this->asImp_().problem().gridGeometry().elementMapper().index(this->element());
106 if (partialReassembler
107 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
108 {
109 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
110
111 for (const auto& localDof : localDofs(this->fvGeometry()))
112 res[localDof.dofIndex()] += residual[localDof.index()];
113
114 // assemble the coupling blocks for coupled models (does nothing if not coupled)
115 maybeAssembleCouplingBlocks(residual);
116 }
117 else if (!this->elementIsGhost())
118 {
119 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
120
121 for (const auto& localDof : localDofs(this->fvGeometry()))
122 res[localDof.dofIndex()] += residual[localDof.index()];
123
124 // assemble the coupling blocks for coupled models (does nothing if not coupled)
125 maybeAssembleCouplingBlocks(residual);
126 }
127 else
128 {
129 // Treatment of ghost elements
130 assert(this->elementIsGhost());
131
132 // handle dofs per codimension
133 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
134 Dune::Hybrid::forEach(std::make_integer_sequence<int, dim+1>{}, [&](auto d)
135 {
136 constexpr int codim = dim - d;
137 const auto& localCoeffs = gridGeometry.feCache().get(this->element().type()).localCoefficients();
138 for (int idx = 0; idx < localCoeffs.size(); ++idx)
139 {
140 const auto& localKey = localCoeffs.localKey(idx);
141
142 // skip if we are not handling this codim right now
143 if (localKey.codim() != codim)
144 continue;
145
146 // do not change the non-ghost entities
147 auto entity = this->element().template subEntity<codim>(localKey.subEntity());
148 if (entity.partitionType() == Dune::InteriorEntity || entity.partitionType() == Dune::BorderEntity)
149 continue;
150
151 // Set identity rows for ALL DOFs of this ghost entity.
152 // Entities with multiple DOFs (e.g. PQ3 edge interior DOFs with 2 per edge)
153 // require iterating over all DOF indices via asMultiMapper(dofMapper()).indices(entity).
154 using BlockType = typename JacobianMatrix::block_type;
155 for (const auto dofIndex : asMultiMapper(gridGeometry.dofMapper()).indices(entity))
156 {
157 BlockType &J = jac[dofIndex][dofIndex];
158 for (int j = 0; j < BlockType::rows; ++j)
159 J[j][j] = 1.0;
160 res[dofIndex] = 0;
161 }
162 }
163 });
164 }
165
166 auto applyDirichlet = [&] (const auto& scvOrLocalDofI,
167 const auto& dirichletValues,
168 const auto eqIdx,
169 const auto pvIdx)
170 {
171 res[scvOrLocalDofI.dofIndex()][eqIdx] = this->curElemVolVars()[scvOrLocalDofI].priVars()[pvIdx] - dirichletValues[pvIdx];
172
173 auto& row = jac[scvOrLocalDofI.dofIndex()];
174 for (auto col = row.begin(); col != row.end(); ++col)
175 row[col.index()][eqIdx] = 0.0;
176
177 jac[scvOrLocalDofI.dofIndex()][scvOrLocalDofI.dofIndex()][eqIdx][pvIdx] = 1.0;
178
179 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
180 if (this->asImp_().problem().gridGeometry().dofOnPeriodicBoundary(scvOrLocalDofI.dofIndex()))
181 {
182 const auto periodicDof = this->asImp_().problem().gridGeometry().periodicallyMappedDof(scvOrLocalDofI.dofIndex());
183 res[periodicDof][eqIdx] = this->asImp_().curSol()[periodicDof][pvIdx] - dirichletValues[pvIdx];
184
185 auto& rowP = jac[periodicDof];
186 for (auto col = rowP.begin(); col != rowP.end(); ++col)
187 rowP[col.index()][eqIdx] = 0.0;
188
189 rowP[periodicDof][eqIdx][pvIdx] = 1.0;
190 }
191 };
192
193 this->asImp_().enforceDirichletConstraints(applyDirichlet);
194 }
195
200 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
201 {
202 this->asImp_().bindLocalViews();
203 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
204
205 auto applyDirichlet = [&] (const auto& scvOrLocalDofI,
206 const auto& dirichletValues,
207 const auto eqIdx,
208 const auto pvIdx)
209 {
210 auto& row = jac[scvOrLocalDofI.dofIndex()];
211 for (auto col = row.begin(); col != row.end(); ++col)
212 row[col.index()][eqIdx] = 0.0;
213
214 jac[scvOrLocalDofI.dofIndex()][scvOrLocalDofI.dofIndex()][eqIdx][pvIdx] = 1.0;
215 };
216
217 this->asImp_().enforceDirichletConstraints(applyDirichlet);
218 }
219
223 template <class ResidualVector>
224 void assembleResidual(ResidualVector& res)
225 {
226 this->asImp_().bindLocalViews();
227 const auto residual = this->evalLocalResidual();
228
229 for (const auto& localDof : localDofs(this->fvGeometry()))
230 res[localDof.dofIndex()] += residual[localDof.index()];
231
232 auto applyDirichlet = [&] (const auto& scvOrLocalDofI,
233 const auto& dirichletValues,
234 const auto eqIdx,
235 const auto pvIdx)
236 {
237 res[scvOrLocalDofI.dofIndex()][eqIdx] = this->curElemVolVars()[scvOrLocalDofI].priVars()[pvIdx] - dirichletValues[pvIdx];
238 };
239
240 this->asImp_().enforceDirichletConstraints(applyDirichlet);
241 }
242
244 template<typename ApplyFunction>
245 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
246 {
247 // enforce Dirichlet boundary conditions
248 this->asImp_().evalDirichletBoundaries(applyDirichlet);
249 // take care of internal Dirichlet constraints (if enabled)
250 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
251 }
252
256 template< typename ApplyDirichletFunctionType >
257 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
258 {
259 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
260 // and set the residual to (privar - dirichletvalue)
261 // when having the new boundary interface Dirichlet conditions are incorporated via constraints
262 if constexpr (!Detail::hasProblemBoundaryTypesForFaceFunction<Problem, FVElementGeometry>())
263 {
264 if (this->elemBcTypes().hasDirichlet())
265 {
266 for (const auto& scvI : scvs(this->fvGeometry()))
267 {
268 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
269 if (bcTypes.hasDirichlet())
270 {
271 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
272
273 // set the Dirichlet conditions in residual and jacobian
274 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
275 {
276 if (bcTypes.isDirichlet(eqIdx))
277 {
278 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
279 assert(0 <= pvIdx && pvIdx < numEq);
280 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
281 }
282 }
283 }
284 }
285 }
286 }
287 }
288
293 template<class... Args>
294 void maybeUpdateCouplingContext(Args&&...) {}
295
300 template<class... Args>
302
303};
304
314template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
316
322template<class TypeTag, class Assembler, class Implementation>
323class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
324: public CVFELocalAssemblerBase<TypeTag, Assembler,
325 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>,
326 true>
327{
333 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
336
337 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
339
340 static constexpr bool enableGridFluxVarsCache
341 = GridVariables::GridFluxVariablesCache::cachingEnabled;
342 static constexpr bool solutionDependentFluxVarsCache
343 = GridVariables::GridFluxVariablesCache::FluxVariablesCache::isSolDependent;
344
345public:
346
348 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
349 using ParentType::ParentType;
350
357 template <class PartialReassembler = DefaultPartialReassembler>
358 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
359 const PartialReassembler* partialReassembler = nullptr)
360 {
361 // get some aliases for convenience
362 const auto& element = this->element();
363 const auto& fvGeometry = this->fvGeometry();
364 const auto& curSol = this->asImp_().curSol();
365
366 auto&& curElemVolVars = this->curElemVolVars();
367 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
368
369 // get the vector of the actual element residuals
370 const auto origResiduals = this->evalLocalResidual();
371
373 // //
374 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
375 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
376 // actual element. In the actual element we evaluate the derivative of the entire residual. //
377 // //
379
380 // if all volvars in the stencil have to be updated or if it's enough to only update the
381 // volVars for the scv whose associated dof has been deflected
382 static const bool updateAllVolVars = getParamFromGroup<bool>(
383 this->asImp_().problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
384 );
385
386 // create the element solution
387 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
388
389 // create the vector storing the partial derivatives
390 ElementResidualVector partialDerivs(Detail::LocalDofs::numLocalDofs(fvGeometry));
391
392 auto deflectionPolicy = Detail::CVFE::makeVariablesDeflectionPolicy(
393 gridVariables.curGridVolVars(),
396 updateAllVolVars
397 );
398
399 auto assembleDerivative = [&, this](const auto& localDof)
400 {
401 // dof index and corresponding actual pri vars
402 const auto dofIdx = localDof.dofIndex();
403 const auto localIdx = localDof.index();
404 deflectionPolicy.store(localDof);
405
406 // calculate derivatives w.r.t to the privars at the dof at hand
407 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
408 {
409 partialDerivs = 0.0;
410
411 auto evalResiduals = [&](PrimaryVariable priVar)
412 {
413 // update the volume variables and compute element residual
414 elemSol[localIdx][pvIdx] = priVar;
415 deflectionPolicy.update(elemSol, localDof, this->asImp_().problem());
416 if constexpr (solutionDependentFluxVarsCache)
417 {
419 if constexpr (enableGridFluxVarsCache)
420 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
421 }
422 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
423 return this->evalLocalResidual();
424 };
425
426 // derive the residuals numerically
427 static const NumericEpsilon<PrimaryVariable, numEq> eps_{this->asImp_().problem().paramGroup()};
428 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
429 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[localIdx][pvIdx], partialDerivs, origResiduals,
430 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
431
432 // update the global stiffness matrix with the current partial derivatives
433 for (const auto& localDofJ : localDofs(fvGeometry))
434 {
435 // don't add derivatives for green dofs
436 if (!partialReassembler
437 || partialReassembler->dofColor(localDofJ.dofIndex()) != EntityColor::green)
438 {
439 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
440 {
441 // A[i][col][eqIdx][pvIdx] is the rate of change of
442 // the residual of equation 'eqIdx' at dof 'i'
443 // depending on the primary variable 'pvIdx' at dof
444 // 'col'.
445 A[localDofJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[localDofJ.index()][eqIdx];
446 }
447 }
448 }
449
450 // restore the original state of the scv's volume variables
451 deflectionPolicy.restore(localDof);
452
453 // restore the original element solution
454 elemSol[localIdx][pvIdx] = curSol[localDof.dofIndex()][pvIdx];
455 this->asImp_().maybeUpdateCouplingContext(localDof, elemSol, pvIdx);
456 }
457 };
458
459 // calculation of the derivatives
460 for (const auto& localDof : localDofs(fvGeometry))
461 assembleDerivative(localDof);
462
463 // restore original state of the flux vars cache in case of global caching.
464 // In the case of local caching this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
465 if constexpr (enableGridFluxVarsCache)
466 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
467
468 // evaluate additional derivatives that might arise from the coupling (no-op if not coupled)
469 this->asImp_().maybeEvalAdditionalDomainDerivatives(origResiduals, A, gridVariables);
470
471 return origResiduals;
472 }
473
474}; // implicit CVFEAssembler with numeric Jacobian
475
481template<class TypeTag, class Assembler, class Implementation>
482class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false, Implementation>
483: public CVFELocalAssemblerBase<TypeTag, Assembler,
484 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false, Implementation>>,
485 false>
486{
493
494 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
496
497public:
498
500 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
501 using ParentType::ParentType;
502
509 template <class PartialReassembler = DefaultPartialReassembler>
510 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
511 const PartialReassembler* partialReassembler = nullptr)
512 {
513 if (partialReassembler)
514 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
515
516 // get some aliases for convenience
517 const auto& element = this->element();
518 const auto& fvGeometry = this->fvGeometry();
519 const auto& curSol = this->asImp_().curSol();
520 auto&& curElemVolVars = this->curElemVolVars();
521
522 // get the vecor of the actual element residuals
523 const auto origResiduals = this->evalLocalResidual();
524 const auto origStorageResiduals = this->evalLocalStorageResidual();
525
527 // //
528 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
529 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
530 // actual element. In the actual element we evaluate the derivative of the entire residual. //
531 // //
533
534 // create the element solution
535 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
536
537 // create the vector storing the partial derivatives
538 ElementResidualVector partialDerivs(Detail::LocalDofs::numLocalDofs(fvGeometry));
539
540 // calculation of the derivatives
541 for (const auto& scv : scvs(fvGeometry))
542 {
543 // dof index and corresponding actual pri vars
544 const auto localIdx = scv.localDofIndex();
545 const auto dofIdx = scv.dofIndex();
546 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
547 const VolumeVariables origVolVars(curVolVars);
548
549 // calculate derivatives w.r.t to the privars at the dof at hand
550 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
551 {
552 partialDerivs = 0.0;
553
554 auto evalStorage = [&](Scalar priVar)
555 {
556 // auto partialDerivsTmp = partialDerivs;
557 elemSol[localIdx][pvIdx] = priVar;
558 curVolVars.update(elemSol, this->asImp_().problem(), element, scv);
559 return this->evalLocalStorageResidual();
560 };
561
562 // derive the residuals numerically
563 static const NumericEpsilon<Scalar, numEq> eps_{this->asImp_().problem().paramGroup()};
564 static const int numDiffMethod = getParamFromGroup<int>(this->asImp_().problem().paramGroup(), "Assembly.NumericDifferenceMethod");
565 NumericDifferentiation::partialDerivative(evalStorage, elemSol[localIdx][pvIdx], partialDerivs, origStorageResiduals,
566 eps_(elemSol[localIdx][pvIdx], pvIdx), numDiffMethod);
567
568 // update the global stiffness matrix with the current partial derivatives
569 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
570 {
571 // A[i][col][eqIdx][pvIdx] is the rate of change of
572 // the residual of equation 'eqIdx' at dof 'i'
573 // depending on the primary variable 'pvIdx' at dof
574 // 'col'.
575 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[localIdx][eqIdx];
576 }
577
578 // restore the original state of the scv's volume variables
579 curVolVars = origVolVars;
580
581 // restore the original element solution
582 elemSol[localIdx][pvIdx] = curSol[dofIdx][pvIdx];
583 // TODO additional dof dependencies
584 }
585 }
586
587 return origResiduals;
588 }
589}; // explicit CVFEAssembler with numeric Jacobian
590
596template<class TypeTag, class Assembler, class Implementation>
597class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true, Implementation>
598: public CVFELocalAssemblerBase<TypeTag, Assembler,
599 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true, Implementation>>,
600 true>
601{
606
607public:
608
610 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
611 using ParentType::ParentType;
612
619 template <class PartialReassembler = DefaultPartialReassembler>
620 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
621 const PartialReassembler* partialReassembler = nullptr)
622 {
623 if (partialReassembler)
624 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
625
626 // get some aliases for convenience
627 const auto& element = this->element();
628 const auto& fvGeometry = this->fvGeometry();
629 const auto& problem = this->asImp_().problem();
630 const auto& curElemVolVars = this->curElemVolVars();
631 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
632
633 // get the vecor of the actual element residuals
634 const auto origResiduals = this->evalLocalResidual();
635
637 // //
638 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
639 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
640 // actual element. In the actual element we evaluate the derivative of the entire residual. //
641 // //
643
644 // calculation of the source and storage derivatives
645 for (const auto& scv : scvs(fvGeometry))
646 {
647 // dof index and corresponding actual pri vars
648 const auto dofIdx = scv.dofIndex();
649 const auto& volVars = curElemVolVars[scv];
650
651 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
652 // only if the problem is instationary we add derivative of storage term
653 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
654 if (!this->assembler().isStationaryProblem())
655 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
656 problem,
657 element,
659 volVars,
660 scv);
661
662 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
663 // add source term derivatives
664 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
665 problem,
666 element,
668 volVars,
669 scv);
670 }
671
672 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
673 for (const auto& scvf : scvfs(fvGeometry))
674 {
675 if (!scvf.boundary())
676 {
677 // add flux term derivatives
678 this->localResidual().addFluxDerivatives(A,
679 problem,
680 element,
684 scvf);
685 }
686
687 // the boundary gets special treatment to simplify
688 // for the user
689 else
690 {
691 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
692 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
693 {
694 // add flux term derivatives
695 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
696 problem,
697 element,
701 scvf);
702 }
703 }
704 }
705
706 return origResiduals;
707 }
708
709}; // implicit CVFEAssembler with analytic Jacobian
710
716template<class TypeTag, class Assembler, class Implementation>
717class CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false, Implementation>
718: public CVFELocalAssemblerBase<TypeTag, Assembler,
719 Detail::CVFE::Impl<Implementation, CVFELocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false, Implementation>>,
720 false>
721{
726
727public:
728
730 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
731 using ParentType::ParentType;
732
739 template <class PartialReassembler = DefaultPartialReassembler>
740 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
741 const PartialReassembler* partialReassembler = nullptr)
742 {
743 if (partialReassembler)
744 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
745
746 // get some aliases for convenience
747 const auto& element = this->element();
748 const auto& fvGeometry = this->fvGeometry();
749 const auto& problem = this->asImp_().problem();
750 const auto& curElemVolVars = this->curElemVolVars();
751
752 // get the vecor of the actual element residuals
753 const auto origResiduals = this->evalLocalResidual();
754
756 // //
757 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
758 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
759 // actual element. In the actual element we evaluate the derivative of the entire residual. //
760 // //
762
763 // calculation of the source and storage derivatives
764 for (const auto& scv : scvs(fvGeometry))
765 {
766 // dof index and corresponding actual pri vars
767 const auto dofIdx = scv.dofIndex();
768 const auto& volVars = curElemVolVars[scv];
769
770 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
771 // only if the problem is instationary we add derivative of storage term
772 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
773 problem,
774 element,
776 volVars,
777 scv);
778 }
779
780 return origResiduals;
781 }
782
783}; // explicit CVFEAssembler with analytic Jacobian
784
785} // end namespace Dumux
786
787#endif
A base class for all local assemblers.
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:730
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:740
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:729
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:620
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:610
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:609
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:348
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:347
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:358
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition assembly/cvfelocalassembler.hh:500
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:510
typename ParentType::LocalResidual LocalResidual
Definition assembly/cvfelocalassembler.hh:499
A base class for all local CVFE assemblers.
Definition assembly/cvfelocalassembler.hh:69
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition assembly/cvfelocalassembler.hh:245
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition assembly/cvfelocalassembler.hh:224
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition assembly/cvfelocalassembler.hh:294
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:200
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition assembly/cvfelocalassembler.hh:257
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:100
void bindLocalViews()
Definition assembly/cvfelocalassembler.hh:88
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition assembly/cvfelocalassembler.hh:301
An assembler for Jacobian and residual contribution per element (CVFE methods).
Definition assembly/cvfelocalassembler.hh:315
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
std::decay_t< decltype(std::declval< Assembler >().localResidual())> LocalResidual
Definition assembly/fvlocalassemblerbase.hh:55
const Element & element() const
Definition assembly/fvlocalassemblerbase.hh:237
VolumeVariables & getVolVarAccess(GridVolumeVariables &gridVolVars, ElementVolumeVariables &elemVolVars, const SubControlVolume &scv)
Definition assembly/fvlocalassemblerbase.hh:304
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:32
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
@ 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
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
Class representing dofs on elements for control-volume finite element schemes.
Adapter to expose a multi-DOF mapper interface for single- and multi-DOF mappers.
Definition elementvariables.hh:24
Definition adapt.hh:17
constexpr auto asMultiMapper(const Mapper &mapper)
Definition multimapperview.hh:48
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition localdof.hh:50
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.