version 3.8
assembly/cclocalassembler.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_CC_LOCAL_ASSEMBLER_HH
14#define DUMUX_CC_LOCAL_ASSEMBLER_HH
15
16#include <dune/common/reservedvector.hh>
17#include <dune/grid/common/gridenums.hh> // for GhostEntity
18#include <dune/istl/matrixindexset.hh>
19
32
33namespace Dumux {
34
44template<class TypeTag, class Assembler, class Implementation, bool implicit>
45class CCLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
46{
51 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
53
54public:
55
56 using ParentType::ParentType;
57
62 template <class ResidualVector, class PartialReassembler = DefaultPartialReassembler>
63 void assembleJacobianAndResidual(JacobianMatrix& jac, ResidualVector& res, GridVariables& gridVariables,
64 const PartialReassembler* partialReassembler)
65 {
66 this->asImp_().bindLocalViews();
67 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
68 if (partialReassembler
69 && partialReassembler->elementColor(globalI) == EntityColor::green)
70 {
71 res[globalI] = this->asImp_().evalLocalResidual()[0]; // forward to the internal implementation
72 }
73 else
74 {
75 res[globalI] = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
76 }
77 }
78
83 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
84 {
85 this->asImp_().bindLocalViews();
86 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
87 }
88
92 template <class ResidualVector>
93 void assembleResidual(ResidualVector& res)
94 {
95 this->asImp_().bindLocalViews();
96 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
97 res[globalI] = this->asImp_().evalLocalResidual()[0]; // forward to the internal implementation
98
100 if constexpr (Problem::enableInternalDirichletConstraints())
101 {
102 const auto applyDirichlet = [&] (const auto& scvI,
103 const auto& dirichletValues,
104 const auto eqIdx,
105 const auto pvIdx)
106 {
107 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
108 };
109
110 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
111 }
112 }
113};
114
123template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true>
125
131template<class TypeTag, class Assembler>
132class CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
133: public CCLocalAssemblerBase<TypeTag, Assembler,
134 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true >
135{
140 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
142 using FVElementGeometry = typename GridGeometry::LocalView;
145 using Problem = typename GridVariables::GridVolumeVariables::Problem;
146
149
151 static constexpr int maxElementStencilSize = GridGeometry::maxElementStencilSize;
153
154public:
155
157
164 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables)
165 {
167 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
168 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
169 // actual element. In the actual element we evaluate the derivative of the entire residual. //
171
172 // get some aliases for convenience
173 const auto& element = this->element();
174 const auto& fvGeometry = this->fvGeometry();
175 const auto& gridGeometry = this->assembler().gridGeometry();
176 auto&& curElemVolVars = this->curElemVolVars();
177 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
178
179 // get stencil information
180 const auto globalI = gridGeometry.elementMapper().index(element);
181 const auto& connectivityMap = gridGeometry.connectivityMap();
182 const auto numNeighbors = connectivityMap[globalI].size();
183
184 // container to store the neighboring elements
185 Dune::ReservedVector<Element, maxElementStencilSize> neighborElements;
186 neighborElements.resize(numNeighbors);
187
188 // assemble the undeflected residual
190 Residuals origResiduals(numNeighbors + 1); origResiduals = 0.0;
191 origResiduals[0] = this->evalLocalResidual()[0];
192
193 // lambda for convenient evaluation of the fluxes across scvfs in the neighbors
194 // if the neighbor is a ghost we don't want to add anything to their residual
195 // so we return 0 and omit computing the flux
196 auto evalNeighborFlux = [&] (const auto& neighbor, const auto& scvf)
197 {
198 if (neighbor.partitionType() == Dune::GhostEntity)
199 return NumEqVector(0.0);
200 else
201 return this->localResidual().evalFlux(this->problem(),
202 neighbor,
203 this->fvGeometry(),
204 this->curElemVolVars(),
205 this->elemFluxVarsCache(), scvf);
206 };
207
208 // get the elements in which we need to evaluate the fluxes
209 // and calculate these in the undeflected state
210 unsigned int j = 1;
211 for (const auto& dataJ : connectivityMap[globalI])
212 {
213 neighborElements[j-1] = gridGeometry.element(dataJ.globalJ);
214 for (const auto scvfIdx : dataJ.scvfsJ)
215 origResiduals[j] += evalNeighborFlux(neighborElements[j-1], fvGeometry.scvf(scvfIdx));
216
217 ++j;
218 }
219
220 // reference to the element's scv (needed later) and corresponding vol vars
221 const auto& scv = fvGeometry.scv(globalI);
222 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
223
224 // save a copy of the original privars and vol vars in order
225 // to restore the original solution after deflection
226 const auto& curSol = this->curSol();
227 const auto origPriVars = curSol[globalI];
228 const auto origVolVars = curVolVars;
229
230 // element solution container to be deflected
231 auto elemSol = elementSolution(element, curSol, gridGeometry);
232
233 // derivatives in the neighbors with respect to the current elements
234 // in index 0 we save the derivative of the element residual with respect to it's own dofs
235 Residuals partialDerivs(numNeighbors + 1);
236
237 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
238 {
239 partialDerivs = 0.0;
240
241 auto evalResiduals = [&](Scalar priVar)
242 {
243 Residuals partialDerivsTmp(numNeighbors + 1);
244 partialDerivsTmp = 0.0;
245 // update the volume variables and the flux var cache
246 elemSol[0][pvIdx] = priVar;
247 curVolVars.update(elemSol, this->problem(), element, scv);
248 elemFluxVarsCache.update(element, fvGeometry, curElemVolVars);
249 if (enableGridFluxVarsCache)
250 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
251
252 // calculate the residual with the deflected primary variables
253 partialDerivsTmp[0] = this->evalLocalResidual()[0];
254
255 // calculate the fluxes in the neighbors with the deflected primary variables
256 for (std::size_t k = 0; k < numNeighbors; ++k)
257 for (auto scvfIdx : connectivityMap[globalI][k].scvfsJ)
258 partialDerivsTmp[k+1] += evalNeighborFlux(neighborElements[k], fvGeometry.scvf(scvfIdx));
259
260 return partialDerivsTmp;
261 };
262
263 // derive the residuals numerically
264 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
265 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
266 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[0][pvIdx], partialDerivs, origResiduals,
267 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
268
269 // Correct derivative for ghost elements, i.e. set a 1 for the derivative w.r.t. the
270 // current primary variable and a 0 elsewhere. As we always solve for a delta of the
271 // solution with respect to the initial one, this results in a delta of 0 for ghosts.
272 if (this->elementIsGhost())
273 {
274 partialDerivs[0] = 0.0;
275 partialDerivs[0][pvIdx] = 1.0;
276 }
277
278 // restore the original state of the scv's volume variables
279 curVolVars = origVolVars;
280
281 // restore the current element solution
282 elemSol[0][pvIdx] = origPriVars[pvIdx];
283
284 // add the current partial derivatives to the global jacobian matrix
285 // no special treatment is needed if globalJ is a ghost because then derivatives have been assembled to 0 above
286 if constexpr (Problem::enableInternalDirichletConstraints())
287 {
288 // check if own element has internal Dirichlet constraint
289 const auto internalDirichletConstraintsOwnElement = this->problem().hasInternalDirichletConstraint(this->element(), scv);
290 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
291
292 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
293 {
294 if (internalDirichletConstraintsOwnElement[eqIdx])
295 {
296 origResiduals[0][eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
297 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
298 }
299 else
300 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
301 }
302
303 // off-diagonal entries
304 j = 1;
305 for (const auto& dataJ : connectivityMap[globalI])
306 {
307 const auto& neighborElement = neighborElements[j-1];
308 const auto& neighborScv = fvGeometry.scv(dataJ.globalJ);
309 const auto internalDirichletConstraintsNeighbor = this->problem().hasInternalDirichletConstraint(neighborElement, neighborScv);
310
311 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
312 {
313 if (internalDirichletConstraintsNeighbor[eqIdx])
314 A[dataJ.globalJ][globalI][eqIdx][pvIdx] = 0.0;
315 else
316 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j][eqIdx];
317 }
318
319 ++j;
320 }
321 }
322 else // no internal Dirichlet constraints specified
323 {
324 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
325 {
326 // the diagonal entries
327 A[globalI][globalI][eqIdx][pvIdx] += partialDerivs[0][eqIdx];
328
329 // off-diagonal entries
330 j = 1;
331 for (const auto& dataJ : connectivityMap[globalI])
332 A[dataJ.globalJ][globalI][eqIdx][pvIdx] += partialDerivs[j++][eqIdx];
333 }
334 }
335 }
336
337 // restore original state of the flux vars cache in case of global caching.
338 // This has to be done in order to guarantee that everything is in an undeflected
339 // state before the assembly of another element is called. In the case of local caching
340 // this is obsolete because the elemFluxVarsCache used here goes out of scope after this.
341 // We only have to do this for the last primary variable, for all others the flux var cache
342 // is updated with the correct element volume variables before residual evaluations
343 if (enableGridFluxVarsCache)
344 gridVariables.gridFluxVarsCache().updateElement(element, fvGeometry, curElemVolVars);
345
346 // return the original residual
347 return origResiduals[0];
348 }
349};
350
351
357template<class TypeTag, class Assembler>
358class CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
359: public CCLocalAssemblerBase<TypeTag, Assembler,
360 CCLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
361{
366 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
369 using Problem = typename GridVariables::GridVolumeVariables::Problem;
370
372
373public:
375
382 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables)
383 {
384 if (this->assembler().isStationaryProblem())
385 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
386
387 // assemble the undeflected residual
388 auto residual = this->evalLocalResidual()[0];
389 const auto storageResidual = this->evalLocalStorageResidual()[0];
390
392 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
393 // neighboring elements all derivatives are zero. For the assembled element only the storage //
394 // derivatives are non-zero. //
396
397 // get some aliases for convenience
398 const auto& element = this->element();
399 const auto& fvGeometry = this->fvGeometry();
400 const auto& gridGeometry = this->assembler().gridGeometry();
401 auto&& curElemVolVars = this->curElemVolVars();
402
403 // reference to the element's scv (needed later) and corresponding vol vars
404 const auto globalI = gridGeometry.elementMapper().index(element);
405 const auto& scv = fvGeometry.scv(globalI);
406 auto& curVolVars = ParentType::getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
407
408 // save a copy of the original privars and vol vars in order
409 // to restore the original solution after deflection
410 const auto& curSol = this->curSol();
411 const auto origPriVars = curSol[globalI];
412 const auto origVolVars = curVolVars;
413
414 // element solution container to be deflected
415 auto elemSol = elementSolution(element, curSol, gridGeometry);
416
417 NumEqVector partialDeriv;
418
419 // derivatives in the neighbors with respect to the current elements
420 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
421 {
422 // reset derivatives of element dof with respect to itself
423 partialDeriv = 0.0;
424
425 auto evalStorage = [&](Scalar priVar)
426 {
427 // update the volume variables and calculate
428 // the residual with the deflected primary variables
429 elemSol[0][pvIdx] = priVar;
430 curVolVars.update(elemSol, this->problem(), element, scv);
431 return this->evalLocalStorageResidual()[0];
432 };
433
434 // for non-ghosts compute the derivative numerically
435 if (!this->elementIsGhost())
436 {
437 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
438 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
439 NumericDifferentiation::partialDerivative(evalStorage, elemSol[0][pvIdx], partialDeriv, storageResidual,
440 eps_(elemSol[0][pvIdx], pvIdx), numDiffMethod);
441 }
442
443 // for ghost elements we assemble a 1.0 where the primary variable and zero everywhere else
444 // as we always solve for a delta of the solution with respect to the initial solution this
445 // results in a delta of zero for ghosts
446 else partialDeriv[pvIdx] = 1.0;
447
448 // restore the original state of the scv's volume variables
449 curVolVars = origVolVars;
450
451 // restore the current element solution
452 elemSol[0][pvIdx] = origPriVars[pvIdx];
453
454 // add the current partial derivatives to the global jacobian matrix
455 // only diagonal entries for explicit jacobians
456 if constexpr (Problem::enableInternalDirichletConstraints())
457 {
458 // check if own element has internal Dirichlet constraint
459 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
460 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
461
462 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
463 {
464 if (internalDirichletConstraints[eqIdx])
465 {
466 residual[eqIdx] = origVolVars.priVars()[eqIdx] - dirichletValues[eqIdx];
467 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
468 }
469 else
470 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
471 }
472 }
473 else
474 {
475 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
476 A[globalI][globalI][eqIdx][pvIdx] += partialDeriv[eqIdx];
477 }
478 }
479
480 // return the original residual
481 return residual;
482 }
483};
484
490template<class TypeTag, class Assembler>
491class CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
492: public CCLocalAssemblerBase<TypeTag, Assembler,
493 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
494{
500 using Problem = typename GridVariables::GridVolumeVariables::Problem;
501
503
504public:
506
513 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, const GridVariables& gridVariables)
514 {
515 // treat ghost separately, we always want zero update for ghosts
516 if (this->elementIsGhost())
517 {
518 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
519 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
520 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
521
522 // return zero residual
523 return NumEqVector(0.0);
524 }
525
526 // assemble the undeflected residual
527 auto residual = this->evalLocalResidual()[0];
528
529 // get some aliases for convenience
530 const auto& problem = this->problem();
531 const auto& element = this->element();
532 const auto& fvGeometry = this->fvGeometry();
533 const auto& curElemVolVars = this->curElemVolVars();
534 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
535
536 // get reference to the element's current vol vars
537 const auto globalI = this->assembler().gridGeometry().elementMapper().index(element);
538 const auto& scv = fvGeometry.scv(globalI);
539 const auto& volVars = curElemVolVars[scv];
540
541 // if the problem is instationary, add derivative of storage term
542 if (!this->assembler().isStationaryProblem())
543 this->localResidual().addStorageDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
544
545 // add source term derivatives
546 this->localResidual().addSourceDerivatives(A[globalI][globalI], problem, element, fvGeometry, volVars, scv);
547
548 // add flux derivatives for each scvf
549 for (const auto& scvf : scvfs(fvGeometry))
550 {
551 // inner faces
552 if (!scvf.boundary())
553 this->localResidual().addFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
554
555 // boundary faces
556 else
557 {
558 const auto& bcTypes = problem.boundaryTypes(element, scvf);
559
560 // add Dirichlet boundary flux derivatives
561 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
562 this->localResidual().addCCDirichletFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
563
564 // add Robin ("solution dependent Neumann") boundary flux derivatives
565 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
566 this->localResidual().addRobinFluxDerivatives(A[globalI], problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
567
568 else
569 DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
570 }
571 }
572
573 if constexpr (Problem::enableInternalDirichletConstraints())
574 {
575 // check if own element has internal Dirichlet constraint
576 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
577 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
578
579 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
580 {
581 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
582 {
583 if (internalDirichletConstraints[eqIdx])
584 {
585 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
586 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
587
588 // inner faces
589 for (const auto& scvf : scvfs(fvGeometry))
590 if (!scvf.boundary())
591 A[globalI][fvGeometry.scv(scvf.outsideScvIdx()).dofIndex()][eqIdx][pvIdx] = 0.0;
592 }
593 }
594 }
595 }
596
597 // return element residual
598 return residual;
599 }
600};
601
607template<class TypeTag, class Assembler>
608class CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
609: public CCLocalAssemblerBase<TypeTag, Assembler,
610 CCLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
611{
617 using Problem = typename GridVariables::GridVolumeVariables::Problem;
618
620
621public:
623
630 NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix& A, const GridVariables& gridVariables)
631 {
632 // treat ghost separately, we always want zero update for ghosts
633 if (this->elementIsGhost())
634 {
635 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
636 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
637 A[globalI][globalI][pvIdx][pvIdx] = 1.0;
638
639 // return zero residual
640 return NumEqVector(0.0);
641 }
642
643 // assemble the undeflected residual
644 const auto residual = this->evalLocalResidual()[0];
645
646 // get reference to the element's current vol vars
647 const auto globalI = this->assembler().gridGeometry().elementMapper().index(this->element());
648 const auto& scv = this->fvGeometry().scv(globalI);
649 const auto& volVars = this->curElemVolVars()[scv];
650
651 // add hand-code derivative of storage term
652 this->localResidual().addStorageDerivatives(A[globalI][globalI], this->problem(), this->element(), this->fvGeometry(), volVars, scv);
653
654 if constexpr (Problem::enableInternalDirichletConstraints())
655 {
656 // check if own element has internal Dirichlet constraint
657 const auto internalDirichletConstraints = this->problem().hasInternalDirichletConstraint(this->element(), scv);
658 const auto dirichletValues = this->problem().internalDirichlet(this->element(), scv);
659
660 for (int pvIdx = 0; pvIdx < numEq; ++pvIdx)
661 {
662 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
663 {
664 if (internalDirichletConstraints[eqIdx])
665 {
666 residual[eqIdx] = volVars.priVars()[eqIdx] - dirichletValues[eqIdx];
667 A[globalI][globalI][eqIdx][pvIdx] = (eqIdx == pvIdx) ? 1.0 : 0.0;
668 }
669 }
670 }
671 }
672
673 // return the original residual
674 return residual;
675 }
676};
677
678} // end namespace Dumux
679
680#endif
The local element solution class for cell-centered methods.
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and explicit time di...
Definition: assembly/cclocalassembler.hh:611
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:630
Cell-centered scheme local assembler using analytic (hand-coded) differentiation and implicit time di...
Definition: assembly/cclocalassembler.hh:494
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, const GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:513
Cell-centered scheme local assembler using numeric differentiation and explicit time discretization.
Definition: assembly/cclocalassembler.hh:361
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:382
Cell-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: assembly/cclocalassembler.hh:135
NumEqVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:164
A base class for all local cell-centered assemblers.
Definition: assembly/cclocalassembler.hh:46
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: assembly/cclocalassembler.hh:83
void assembleResidual(ResidualVector &res)
Assemble the residual only.
Definition: assembly/cclocalassembler.hh:93
void assembleJacobianAndResidual(JacobianMatrix &jac, ResidualVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: assembly/cclocalassembler.hh:63
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: assembly/cclocalassembler.hh:124
A base class for all local assemblers.
Definition: assembly/fvlocalassemblerbase.hh:36
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: assembly/fvlocalassemblerbase.hh:253
Implementation & asImp_()
Definition: assembly/fvlocalassemblerbase.hh:297
const Assembler & assembler() const
The assembler.
Definition: assembly/fvlocalassemblerbase.hh:233
const Element & element() const
The current element.
Definition: assembly/fvlocalassemblerbase.hh:237
The flux stencil specialized for different discretization schemes.
Definition: fluxstencil.hh:33
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
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:491
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
Defines all properties used in Dumux.
An enum class to define various differentiation methods available in order to compute the derivatives...
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
The flux stencil specialized for different discretization schemes.
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
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
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 helper to deduce a vector with the same size as numbers of equations.
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.
A arithmetic block vector type based on DUNE's reserved vector.