3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boxlocalassembler.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_BOX_LOCAL_ASSEMBLER_HH
26#define DUMUX_BOX_LOCAL_ASSEMBLER_HH
27
28#include <dune/istl/matrixindexset.hh>
29#include <dune/istl/bvector.hh>
30
40
41namespace Dumux {
42
52template<class TypeTag, class Assembler, class Implementation, bool implicit>
53class BoxLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
54{
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60
62
63public:
64
65 using ParentType::ParentType;
66
68 {
70 this->elemBcTypes().update(this->problem(), this->element(), this->fvGeometry());
71 }
72
73
78 template <class PartialReassembler = DefaultPartialReassembler>
79 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
80 const PartialReassembler* partialReassembler = nullptr)
81 {
82 this->asImp_().bindLocalViews();
83 const auto eIdxGlobal = this->assembler().gridGeometry().elementMapper().index(this->element());
84 if (partialReassembler
85 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
86 {
87 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
88 for (const auto& scv : scvs(this->fvGeometry()))
89 res[scv.dofIndex()] += residual[scv.localDofIndex()];
90 }
91 else if (!this->elementIsGhost())
92 {
93 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
94 for (const auto& scv : scvs(this->fvGeometry()))
95 res[scv.dofIndex()] += residual[scv.localDofIndex()];
96 }
97 else
98 {
99 using GridGeometry = typename GridVariables::GridGeometry;
100 using GridView = typename GridGeometry::GridView;
101 static constexpr auto dim = GridView::dimension;
102
103 int numVerticesLocal = this->element().subEntities(dim);
104
105 for (int i = 0; i < numVerticesLocal; ++i) {
106 auto vertex = this->element().template subEntity<dim>(i);
107
108 if (vertex.partitionType() == Dune::InteriorEntity ||
109 vertex.partitionType() == Dune::BorderEntity)
110 {
111 // do not change the non-ghost vertices
112 continue;
113 }
114
115 // set main diagonal entries for the vertex
116 int vIdx = this->assembler().gridGeometry().vertexMapper().index(vertex);
117
118 typedef typename JacobianMatrix::block_type BlockType;
119 BlockType &J = jac[vIdx][vIdx];
120 for (int j = 0; j < BlockType::rows; ++j)
121 J[j][j] = 1.0;
122
123 // set residual for the vertex
124 res[vIdx] = 0;
125 }
126 }
127
128 auto applyDirichlet = [&] (const auto& scvI,
129 const auto& dirichletValues,
130 const auto eqIdx,
131 const auto pvIdx)
132 {
133 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
134
135 auto& row = jac[scvI.dofIndex()];
136 for (auto col = row.begin(); col != row.end(); ++col)
137 row[col.index()][eqIdx] = 0.0;
138
139 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
140
141 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
142 if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
143 {
144 const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
145 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
146 const auto end = jac[periodicDof].end();
147 for (auto it = jac[periodicDof].begin(); it != end; ++it)
148 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
149 }
150 };
151
152 this->asImp_().enforceDirichletConstraints(applyDirichlet);
153 }
154
159 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
160 {
161 this->asImp_().bindLocalViews();
162 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
163
164 auto applyDirichlet = [&] (const auto& scvI,
165 const auto& dirichletValues,
166 const auto eqIdx,
167 const auto pvIdx)
168 {
169 auto& row = jac[scvI.dofIndex()];
170 for (auto col = row.begin(); col != row.end(); ++col)
171 row[col.index()][eqIdx] = 0.0;
172
173 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
174 };
175
176 this->asImp_().enforceDirichletConstraints(applyDirichlet);
177 }
178
182 void assembleResidual(SolutionVector& res)
183 {
184 this->asImp_().bindLocalViews();
185 const auto residual = this->evalLocalResidual();
186
187 for (const auto& scv : scvs(this->fvGeometry()))
188 res[scv.dofIndex()] += residual[scv.localDofIndex()];
189
190 auto applyDirichlet = [&] (const auto& scvI,
191 const auto& dirichletValues,
192 const auto eqIdx,
193 const auto pvIdx)
194 {
195 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
196 };
197
198 this->asImp_().enforceDirichletConstraints(applyDirichlet);
199 }
200
202 template<typename ApplyFunction>
203 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
204 {
205 // enforce Dirichlet boundary conditions
206 this->asImp_().evalDirichletBoundaries(applyDirichlet);
207 // take care of internal Dirichlet constraints (if enabled)
208 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
209 }
210
214 template< typename ApplyDirichletFunctionType >
215 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
216 {
217 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
218 // and set the residual to (privar - dirichletvalue)
219 if (this->elemBcTypes().hasDirichlet())
220 {
221 for (const auto& scvI : scvs(this->fvGeometry()))
222 {
223 const auto bcTypes = this->elemBcTypes()[scvI.localDofIndex()];
224 if (bcTypes.hasDirichlet())
225 {
226 const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
227
228 // set the Dirichlet conditions in residual and jacobian
229 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
230 {
231 if (bcTypes.isDirichlet(eqIdx))
232 {
233 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
234 assert(0 <= pvIdx && pvIdx < numEq);
235 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
236 }
237 }
238 }
239 }
240 }
241 }
242
243};
244
253template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true>
255
261template<class TypeTag, class Assembler>
262class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
263: public BoxLocalAssemblerBase<TypeTag, Assembler,
264 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true>
265{
273 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
274
277
278 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
279
280public:
281
283
290 template <class PartialReassembler = DefaultPartialReassembler>
291 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
292 const PartialReassembler* partialReassembler = nullptr)
293 {
294 // get some aliases for convenience
295 const auto& element = this->element();
296 const auto& fvGeometry = this->fvGeometry();
297 const auto& curSol = this->curSol();
298 auto&& curElemVolVars = this->curElemVolVars();
299
300 // get the vector of the actual element residuals
301 const auto origResiduals = this->evalLocalResidual();
302
304 // //
305 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
306 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
307 // actual element. In the actual element we evaluate the derivative of the entire residual. //
308 // //
310
311 // create the element solution
312 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
313
314 // create the vector storing the partial derivatives
315 ElementResidualVector partialDerivs(element.subEntities(dim));
316
317 // calculation of the derivatives
318 for (auto&& scv : scvs(fvGeometry))
319 {
320 // dof index and corresponding actual pri vars
321 const auto dofIdx = scv.dofIndex();
322 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
323 const VolumeVariables origVolVars(curVolVars);
324
325 // calculate derivatives w.r.t to the privars at the dof at hand
326 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
327 {
328 partialDerivs = 0.0;
329
330 auto evalResiduals = [&](Scalar priVar)
331 {
332 // update the volume variables and compute element residual
333 elemSol[scv.localDofIndex()][pvIdx] = priVar;
334 curVolVars.update(elemSol, this->problem(), element, scv);
335 return this->evalLocalResidual();
336 };
337
338 // derive the residuals numerically
339 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
340 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
341 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
342 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
343
344 // update the global stiffness matrix with the current partial derivatives
345 for (auto&& scvJ : scvs(fvGeometry))
346 {
347 // don't add derivatives for green vertices
348 if (!partialReassembler
349 || partialReassembler->vertexColor(scvJ.dofIndex()) != EntityColor::green)
350 {
351 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
352 {
353 // A[i][col][eqIdx][pvIdx] is the rate of change of
354 // the residual of equation 'eqIdx' at dof 'i'
355 // depending on the primary variable 'pvIdx' at dof
356 // 'col'.
357 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
358 }
359 }
360 }
361
362 // restore the original state of the scv's volume variables
363 curVolVars = origVolVars;
364
365 // restore the original element solution
366 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
367 // TODO additional dof dependencies
368 }
369 }
370 return origResiduals;
371 }
372
373}; // implicit BoxAssembler with numeric Jacobian
374
380template<class TypeTag, class Assembler>
381class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
382: public BoxLocalAssemblerBase<TypeTag, Assembler,
383 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
384{
392 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
393
396
397public:
398
400
407 template <class PartialReassembler = DefaultPartialReassembler>
408 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
409 const PartialReassembler* partialReassembler = nullptr)
410 {
411 if (partialReassembler)
412 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
413
414 // get some aliases for convenience
415 const auto& element = this->element();
416 const auto& fvGeometry = this->fvGeometry();
417 const auto& curSol = this->curSol();
418 auto&& curElemVolVars = this->curElemVolVars();
419
420 // get the vecor of the acutal element residuals
421 const auto origResiduals = this->evalLocalResidual();
422 const auto origStorageResiduals = this->evalLocalStorageResidual();
423
425 // //
426 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
427 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
428 // actual element. In the actual element we evaluate the derivative of the entire residual. //
429 // //
431
432 // create the element solution
433 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
434
435 // create the vector storing the partial derivatives
436 ElementResidualVector partialDerivs(element.subEntities(dim));
437
438 // calculation of the derivatives
439 for (auto&& scv : scvs(fvGeometry))
440 {
441 // dof index and corresponding actual pri vars
442 const auto dofIdx = scv.dofIndex();
443 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
444 const VolumeVariables origVolVars(curVolVars);
445
446 // calculate derivatives w.r.t to the privars at the dof at hand
447 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
448 {
449 partialDerivs = 0.0;
450
451 auto evalStorage = [&](Scalar priVar)
452 {
453 // auto partialDerivsTmp = partialDerivs;
454 elemSol[scv.localDofIndex()][pvIdx] = priVar;
455 curVolVars.update(elemSol, this->problem(), element, scv);
456 return this->evalLocalStorageResidual();
457 };
458
459 // derive the residuals numerically
460 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
461 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
462 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
463 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
464
465 // update the global stiffness matrix with the current partial derivatives
466 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
467 {
468 // A[i][col][eqIdx][pvIdx] is the rate of change of
469 // the residual of equation 'eqIdx' at dof 'i'
470 // depending on the primary variable 'pvIdx' at dof
471 // 'col'.
472 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
473 }
474
475 // restore the original state of the scv's volume variables
476 curVolVars = origVolVars;
477
478 // restore the original element solution
479 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
480 // TODO additional dof dependencies
481 }
482 }
483 return origResiduals;
484 }
485}; // explicit BoxAssembler with numeric Jacobian
486
492template<class TypeTag, class Assembler>
493class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
494: public BoxLocalAssemblerBase<TypeTag, Assembler,
495 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
496{
502 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
503
504public:
505
507
514 template <class PartialReassembler = DefaultPartialReassembler>
515 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
516 const PartialReassembler* partialReassembler = nullptr)
517 {
518 if (partialReassembler)
519 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
520
521 // get some aliases for convenience
522 const auto& element = this->element();
523 const auto& fvGeometry = this->fvGeometry();
524 const auto& problem = this->problem();
525 const auto& curElemVolVars = this->curElemVolVars();
526 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
527
528 // get the vecor of the acutal element residuals
529 const auto origResiduals = this->evalLocalResidual();
530
532 // //
533 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
534 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
535 // actual element. In the actual element we evaluate the derivative of the entire residual. //
536 // //
538
539 // calculation of the source and storage derivatives
540 for (const auto& scv : scvs(fvGeometry))
541 {
542 // dof index and corresponding actual pri vars
543 const auto dofIdx = scv.dofIndex();
544 const auto& volVars = curElemVolVars[scv];
545
546 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
547 // only if the problem is instationary we add derivative of storage term
548 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
549 if (!this->assembler().isStationaryProblem())
550 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
551 problem,
552 element,
553 fvGeometry,
554 volVars,
555 scv);
556
557 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
558 // add source term derivatives
559 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
560 problem,
561 element,
562 fvGeometry,
563 volVars,
564 scv);
565 }
566
567 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
568 for (const auto& scvf : scvfs(fvGeometry))
569 {
570 if (!scvf.boundary())
571 {
572 // add flux term derivatives
573 this->localResidual().addFluxDerivatives(A,
574 problem,
575 element,
576 fvGeometry,
577 curElemVolVars,
578 elemFluxVarsCache,
579 scvf);
580 }
581
582 // the boundary gets special treatment to simplify
583 // for the user
584 else
585 {
586 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
587 if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
588 {
589 // add flux term derivatives
590 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
591 problem,
592 element,
593 fvGeometry,
594 curElemVolVars,
595 elemFluxVarsCache,
596 scvf);
597 }
598 }
599 }
600
601 return origResiduals;
602 }
603
604}; // implicit BoxAssembler with analytic Jacobian
605
611template<class TypeTag, class Assembler>
612class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
613: public BoxLocalAssemblerBase<TypeTag, Assembler,
614 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
615{
621 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
622
623public:
624
626
633 template <class PartialReassembler = DefaultPartialReassembler>
634 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
635 const PartialReassembler* partialReassembler = nullptr)
636 {
637 if (partialReassembler)
638 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
639
640 // get some aliases for convenience
641 const auto& element = this->element();
642 const auto& fvGeometry = this->fvGeometry();
643 const auto& problem = this->problem();
644 const auto& curElemVolVars = this->curElemVolVars();
645
646 // get the vecor of the acutal element residuals
647 const auto origResiduals = this->evalLocalResidual();
648
650 // //
651 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
652 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
653 // actual element. In the actual element we evaluate the derivative of the entire residual. //
654 // //
656
657 // calculation of the source and storage derivatives
658 for (const auto& scv : scvs(fvGeometry))
659 {
660 // dof index and corresponding actual pri vars
661 const auto dofIdx = scv.dofIndex();
662 const auto& volVars = curElemVolVars[scv];
663
664 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
665 // only if the problem is instationary we add derivative of storage term
666 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
667 problem,
668 element,
669 fvGeometry,
670 volVars,
671 scv);
672 }
673
674 return origResiduals;
675 }
676
677}; // explicit BoxAssembler with analytic Jacobian
678
679} // end namespace Dumux
680
681#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.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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 box assemblers.
Definition: boxlocalassembler.hh:54
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: boxlocalassembler.hh:159
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: boxlocalassembler.hh:215
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: boxlocalassembler.hh:203
void bindLocalViews()
Definition: boxlocalassembler.hh:67
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: boxlocalassembler.hh:182
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: boxlocalassembler.hh:79
An assembler for Jacobian and residual contribution per element (box methods)
Definition: boxlocalassembler.hh:254
Box local assembler using numeric differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:265
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: boxlocalassembler.hh:291
Box local assembler using numeric differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:384
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: boxlocalassembler.hh:408
Box local assembler using analytic differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:496
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: boxlocalassembler.hh:515
Box local assembler using analytic differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:615
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: boxlocalassembler.hh:634
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:263
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:279
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:307
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:239
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:259
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:243
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:251
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:247
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:424
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
Declares all properties used in Dumux.
The local element solution class for the box method.