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