3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
fcdiamondlocalassembler.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_FC_DIAMOND_LOCAL_ASSEMBLER_HH
26#define DUMUX_FC_DIAMOND_LOCAL_ASSEMBLER_HH
27
28#include <dune/grid/common/gridenums.hh>
29
39
40#include "volvardeflectionhelper_.hh"
41
42namespace Dumux {
43
44namespace Detail {
45
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
66template<class TypeTag, class Assembler, class Implementation, bool implicit>
67class FaceCenteredDiamondLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
68{
75
76 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
77
78public:
79
80 using ParentType::ParentType;
81
83 {
85 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
86 }
87
92 template <class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::NoOperator>
93 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
94 const PartialReassembler* partialReassembler,
95 const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction())
96 {
97 static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(), "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!");
98
99 this->asImp_().bindLocalViews();
100 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
101 const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element());
102 if (partialReassembler
103 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
104 {
105 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
106 for (const auto& scv : scvs(this->fvGeometry()))
107 res[scv.dofIndex()] += residual[scv.localDofIndex()];
108
109 // assemble the coupling blocks for coupled models (does nothing if not coupled)
110 maybeAssembleCouplingBlocks(residual);
111 }
112 else if (!this->elementIsGhost())
113 {
114 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
115 for (const auto& scv : scvs(this->fvGeometry()))
116 res[scv.dofIndex()] += residual[scv.localDofIndex()];
117
118 // assemble the coupling blocks for coupled models (does nothing if not coupled)
119 maybeAssembleCouplingBlocks(residual);
120 }
121 else
122 {
123 const auto numLocalFaces = this->element().subEntities(1);
124 for (int i = 0; i < numLocalFaces; ++i)
125 {
126 // do not change the non-ghost entities
127 const auto face = this->element().template subEntity<1>(i);
128 if (face.partitionType() == Dune::InteriorEntity || face.partitionType() == Dune::BorderEntity)
129 continue;
130
131 // set main diagonal entries for the face
132 const auto dofIndex = gridGeometry.dofMapper().index(face);
133
134 typedef typename JacobianMatrix::block_type BlockType;
135 BlockType &J = jac[dofIndex][dofIndex];
136 for (int j = 0; j < BlockType::rows; ++j)
137 J[j][j] = 1.0;
138
139 // set residual for the face
140 res[dofIndex] = 0;
141 }
142 }
143
144 auto applyDirichlet = [&] (const auto& scvI,
145 const auto& dirichletValues,
146 const auto eqIdx,
147 const auto pvIdx)
148 {
149 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
150
151 auto& row = jac[scvI.dofIndex()];
152 for (auto col = row.begin(); col != row.end(); ++col)
153 row[col.index()][eqIdx] = 0.0;
154
155 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
156
157
158 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof TODO periodic
159 // if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
160 // {
161 // const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
162 // res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
163 // const auto end = jac[periodicDof].end();
164 // for (auto it = jac[periodicDof].begin(); it != end; ++it)
165 // (*it) = periodicDof != it.index() ? 0.0 : 1.0;
166 // }
167 };
168
169 this->asImp_().enforceDirichletConstraints(applyDirichlet);
170 }
171
176 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
177 {
178 this->asImp_().bindLocalViews();
179 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
180
181 auto applyDirichlet = [&] (const auto& scvI,
182 const auto& dirichletValues,
183 const auto eqIdx,
184 const auto pvIdx)
185 {
186 auto& row = jac[scvI.dofIndex()];
187 for (auto col = row.begin(); col != row.end(); ++col)
188 row[col.index()][eqIdx] = 0.0;
189
190 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
191 };
192
193 this->asImp_().enforceDirichletConstraints(applyDirichlet);
194 }
195
199 void assembleResidual(SolutionVector& res)
200 {
201 this->asImp_().bindLocalViews();
202 const auto residual = this->evalLocalResidual();
203
204 for (const auto& scv : scvs(this->fvGeometry()))
205 res[scv.dofIndex()] += residual[scv.localDofIndex()];
206
207 auto applyDirichlet = [&] (const auto& scvI,
208 const auto& dirichletValues,
209 const auto eqIdx,
210 const auto pvIdx)
211 {
212 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
213 };
214
215 this->asImp_().enforceDirichletConstraints(applyDirichlet);
216 }
217
221 template<typename ApplyFunction>
222 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
223 {
224 // enforce Dirichlet boundary conditions
225 this->asImp_().evalDirichletBoundaries(applyDirichlet);
226 // take care of internal Dirichlet constraints (if enabled)
227 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
228 }
229
233 template< typename ApplyDirichletFunctionType >
234 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
235 {
236 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
237 // and set the residual to (privar - dirichletvalue)
238 if (this->elemBcTypes().hasDirichlet())
239 {
240 for (const auto& scvf : scvfs(this->fvGeometry()))
241 {
242 if (scvf.boundary())
243 {
244 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvf);
245 if (bcTypes.hasDirichlet())
246 {
247 const PrimaryVariables dirichletValues(this->asImp_().problem().dirichlet(this->element(), scvf));
248 const auto& scv = this->fvGeometry().scv(scvf.insideScvIdx());
249
250 // set the Dirichlet conditions in residual and jacobian
251 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
252 {
253 if (bcTypes.isDirichlet(eqIdx))
254 {
255 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
256 assert(0 <= pvIdx && pvIdx < numEq);
257 applyDirichlet(scv, dirichletValues, eqIdx, pvIdx);
258 }
259 }
260 }
261 }
262 }
263 }
264 }
265
270 template<class... Args>
271 void maybeUpdateCouplingContext(Args&&...) {}
272
277 template<class... Args>
279};
280
289template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
291
297template<class TypeTag, class Assembler, class Implementation>
298class FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
299: public FaceCenteredDiamondLocalAssemblerBase<TypeTag, Assembler,
300 Detail::Impl<Implementation, FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>, true>
301{
305 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
307 using FVElementGeometry = typename GridGeometry::LocalView;
308 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
309 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
313
314 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
316
317public:
318
320 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
322
329 template <class PartialReassembler = DefaultPartialReassembler>
330 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
331 const PartialReassembler* partialReassembler = nullptr)
332 {
333 // get some aliases for convenience
334 const auto& problem = this->asImp_().problem();
335 const auto& element = this->element();
336 const auto& fvGeometry = this->fvGeometry();
337 const auto& curSol = this->asImp_().curSol();
338
339 auto&& curElemVolVars = this->curElemVolVars();
340
341 // get the vector of the actual element residuals
342 const auto origResiduals = this->evalLocalResidual();
343
345 // Calculate derivatives of all dofs in the element with respect to the dofs in the stencil. //
347
348 // if all volvars in the stencil have to be updated or if it's enough to only update the
349 // volVars for the scv whose associated dof has been deflected
350 static const bool updateAllVolVars = getParamFromGroup<bool>(
351 problem.paramGroup(), "Assembly.FCDiamondVolVarsDependOnAllElementDofs", false
352 );
353
354 // create the element solution
355 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
356
357 // one residual per element facet
358 const auto numElementResiduals = fvGeometry.numScv();
359
360 // create the vector storing the partial derivatives
361 ElementResidualVector partialDerivs(numElementResiduals);
362
363 Detail::VolVarsDeflectionHelper deflectionHelper(
364 [&] (const auto& scv) -> VolumeVariables& {
365 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
366 },
367 fvGeometry,
368 updateAllVolVars
369 );
370
371 // calculation of the derivatives
372 for (const auto& scv : scvs(fvGeometry))
373 {
374 // dof index and corresponding actual pri vars
375 const auto dofIdx = scv.dofIndex();
376 deflectionHelper.setCurrent(scv);
377
378 // calculate derivatives w.r.t to the privars at the dof at hand
379 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
380 {
381 partialDerivs = 0.0;
382
383 auto evalResiduals = [&](Scalar priVar)
384 {
385 // update the volume variables and compute element residual
386 elemSol[scv.localDofIndex()][pvIdx] = priVar;
387 deflectionHelper.deflect(elemSol, scv, problem);
388 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
389 return this->evalLocalResidual();
390 };
391
392 // derive the residuals numerically
393 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
394 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
395 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
396 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
397
398 // update the global stiffness matrix with the current partial derivatives
399 for (const auto& scvJ : scvs(fvGeometry))
400 {
401 // TODO partial reassembly
402 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
403 {
404 // A[i][col][eqIdx][pvIdx] is the rate of change of
405 // the residual of equation 'eqIdx' at dof 'i'
406 // depending on the primary variable 'pvIdx' at dof
407 // 'col'.
408 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
409 }
410 }
411
412 // restore the original state of the scv's volume variables
413 deflectionHelper.restore(scv);
414
415 // restore the original element solution
416 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
417 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
418 }
419 }
420
421 return origResiduals;
422 }
423};
424
425
431template<class TypeTag, class Assembler>
432class FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
433: public FaceCenteredDiamondLocalAssemblerBase<TypeTag, Assembler,
434 FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
435{
439 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
443 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
445
447
448public:
450
457 template <class PartialReassembler = DefaultPartialReassembler>
458 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
459 const PartialReassembler* partialReassembler = nullptr)
460 {
461 if (partialReassembler)
462 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
463
464 // get some aliases for convenience
465 const auto& problem = this->asImp_().problem();
466 const auto& element = this->element();
467 const auto& fvGeometry = this->fvGeometry();
468 const auto& curSol = this->asImp_().curSol();
469 auto&& curElemVolVars = this->curElemVolVars();
470
471 // get the vecor of the actual element residuals
472 const auto origResiduals = this->evalLocalResidual();
473 const auto origStorageResiduals = this->evalLocalStorageResidual();
474
476 // //
477 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
478 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
479 // actual element. In the actual element we evaluate the derivative of the entire residual. //
480 // //
482
483 // create the element solution
484 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
485
486 // create the vector storing the partial derivatives
487 ElementResidualVector partialDerivs(fvGeometry.numScv());
488
489 // calculation of the derivatives
490 for (auto&& scv : scvs(fvGeometry))
491 {
492 // dof index and corresponding actual pri vars
493 const auto dofIdx = scv.dofIndex();
494 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
495 const VolumeVariables origVolVars(curVolVars);
496
497 // calculate derivatives w.r.t to the privars at the dof at hand
498 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
499 {
500 partialDerivs = 0.0;
501
502 auto evalStorage = [&](Scalar priVar)
503 {
504 // auto partialDerivsTmp = partialDerivs;
505 elemSol[scv.localDofIndex()][pvIdx] = priVar;
506 curVolVars.update(elemSol, problem, element, scv);
507 return this->evalLocalStorageResidual();
508 };
509
510 // derive the residuals numerically
511 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
512 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
513 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
514 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
515
516 // update the global stiffness matrix with the current partial derivatives
517 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
518 {
519 // A[i][col][eqIdx][pvIdx] is the rate of change of
520 // the residual of equation 'eqIdx' at dof 'i'
521 // depending on the primary variable 'pvIdx' at dof
522 // 'col'.
523 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
524 }
525
526 // restore the original state of the scv's volume variables
527 curVolVars = origVolVars;
528
529 // restore the original element solution
530 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
531 // TODO additional dof dependencies
532 }
533 }
534 return origResiduals;
535 }
536};
537
543template<class TypeTag, class Assembler>
544class FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
545: public FaceCenteredDiamondLocalAssemblerBase<TypeTag, Assembler,
546 FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
547{
553 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
555
557
558public:
560
567 template <class PartialReassembler = DefaultPartialReassembler>
568 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
569 const PartialReassembler* partialReassembler = nullptr)
570 {
571 if (partialReassembler)
572 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
573
574 // get some aliases for convenience
575 const auto& element = this->element();
576 const auto& fvGeometry = this->fvGeometry();
577 const auto& problem = this->asImp_().problem();
578 const auto& curElemVolVars = this->curElemVolVars();
579 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
580
581 // get the vecor of the actual element residuals
582 const auto origResiduals = this->evalLocalResidual();
583
585 // //
586 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
587 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
588 // actual element. In the actual element we evaluate the derivative of the entire residual. //
589 // //
591
592 // calculation of the source and storage derivatives
593 for (const auto& scv : scvs(fvGeometry))
594 {
595 // dof index and corresponding actual pri vars
596 const auto dofIdx = scv.dofIndex();
597 const auto& volVars = curElemVolVars[scv];
598
599 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
600 // only if the problem is instationary we add derivative of storage term
601 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
602 if (!this->assembler().isStationaryProblem())
603 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
604 problem,
605 element,
606 fvGeometry,
607 volVars,
608 scv);
609
610 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
611 // add source term derivatives
612 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
613 problem,
614 element,
615 fvGeometry,
616 volVars,
617 scv);
618 }
619
620 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
621 for (const auto& scvf : scvfs(fvGeometry))
622 {
623 if (!scvf.boundary())
624 {
625 // add flux term derivatives
626 this->localResidual().addFluxDerivatives(A,
627 problem,
628 element,
629 fvGeometry,
630 curElemVolVars,
631 elemFluxVarsCache,
632 scvf);
633 }
634
635 // the boundary gets special treatment to simplify
636 // for the user
637 else
638 {
639 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
640 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
641 {
642 // add flux term derivatives
643 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
644 problem,
645 element,
646 fvGeometry,
647 curElemVolVars,
648 elemFluxVarsCache,
649 scvf);
650 }
651 }
652 }
653
654 return origResiduals;
655 }
656};
657
663template<class TypeTag, class Assembler>
664class FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
665: public FaceCenteredDiamondLocalAssemblerBase<TypeTag, Assembler,
666 FaceCenteredDiamondLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
667{
673 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
675
677
678public:
680
687 template <class PartialReassembler = DefaultPartialReassembler>
688 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
689 const PartialReassembler* partialReassembler = nullptr)
690 {
691 if (partialReassembler)
692 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
693
694 // get some aliases for convenience
695 const auto& element = this->element();
696 const auto& fvGeometry = this->fvGeometry();
697 const auto& problem = this->asImp_().problem();
698 const auto& curElemVolVars = this->curElemVolVars();
699
700 // get the vector of the actual element residuals
701 const auto origResiduals = this->evalLocalResidual();
702
704 // //
705 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
706 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
707 // actual element. In the actual element we evaluate the derivative of the entire residual. //
708 // //
710
711 // calculation of the source and storage derivatives
712 for (const auto& scv : scvs(fvGeometry))
713 {
714 // dof index and corresponding actual pri vars
715 const auto dofIdx = scv.dofIndex();
716 const auto& volVars = curElemVolVars[scv];
717
718 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
719 // only if the problem is instationary we add derivative of storage term
720 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
721 problem,
722 element,
723 fvGeometry,
724 volVars,
725 scv);
726 }
727
728 return origResiduals;
729 }
730};
731
732} // end namespace Dumux
733
734#endif
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.
An adapter class for local assemblers using numeric differentiation.
Detects which entries in the Jacobian have to be recomputed.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
@ green
does not need to be reassembled
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
Definition: fcdiamondlocalassembler.hh:47
constexpr void operator()(Args &&...) const
Definition: fcdiamondlocalassembler.hh:49
A base class for all local face-centered assemblers.
Definition: fcdiamondlocalassembler.hh:68
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: fcdiamondlocalassembler.hh:278
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: fcdiamondlocalassembler.hh:271
void bindLocalViews()
Definition: fcdiamondlocalassembler.hh:82
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: fcdiamondlocalassembler.hh:222
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: fcdiamondlocalassembler.hh:234
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: fcdiamondlocalassembler.hh:176
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler, const CouplingFunction &maybeAssembleCouplingBlocks=CouplingFunction())
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: fcdiamondlocalassembler.hh:93
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: fcdiamondlocalassembler.hh:199
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fcdiamondlocalassembler.hh:290
Face-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: fcdiamondlocalassembler.hh:301
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: fcdiamondlocalassembler.hh:330
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: fcdiamondlocalassembler.hh:319
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fcdiamondlocalassembler.hh:320
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: fcdiamondlocalassembler.hh:458
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: fcdiamondlocalassembler.hh:568
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: fcdiamondlocalassembler.hh:688
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: fvlocalassemblerbase.hh:185
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
const Problem & problem() const
The problem.
Definition: fvlocalassemblerbase.hh:241
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:503
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:61
Declares all properties used in Dumux.