3.6-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/common/reservedvector.hh>
29#include <dune/istl/matrixindexset.hh>
30#include <dune/istl/bvector.hh>
31
41
42#include "volvardeflectionhelper_.hh"
43
44namespace Dumux {
45
55template<class TypeTag, class Assembler, class Implementation, bool implicit>
56class BoxLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
57{
62 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
63
65
66public:
67
68 using ParentType::ParentType;
69
71 {
73 this->elemBcTypes().update(this->problem(), this->element(), this->fvGeometry());
74 }
75
76
81 template <class PartialReassembler = DefaultPartialReassembler>
82 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
83 const PartialReassembler* partialReassembler = nullptr)
84 {
85 this->asImp_().bindLocalViews();
86 const auto eIdxGlobal = this->assembler().gridGeometry().elementMapper().index(this->element());
87 if (partialReassembler
88 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
89 {
90 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
91 for (const auto& scv : scvs(this->fvGeometry()))
92 res[scv.dofIndex()] += residual[scv.localDofIndex()];
93 }
94 else if (!this->elementIsGhost())
95 {
96 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
97 for (const auto& scv : scvs(this->fvGeometry()))
98 res[scv.dofIndex()] += residual[scv.localDofIndex()];
99 }
100 else
101 {
102 using GridGeometry = typename GridVariables::GridGeometry;
103 using GridView = typename GridGeometry::GridView;
104 static constexpr auto dim = GridView::dimension;
105
106 int numVerticesLocal = this->element().subEntities(dim);
107
108 for (int i = 0; i < numVerticesLocal; ++i) {
109 auto vertex = this->element().template subEntity<dim>(i);
110
111 if (vertex.partitionType() == Dune::InteriorEntity ||
112 vertex.partitionType() == Dune::BorderEntity)
113 {
114 // do not change the non-ghost vertices
115 continue;
116 }
117
118 // set main diagonal entries for the vertex
119 int vIdx = this->assembler().gridGeometry().vertexMapper().index(vertex);
120
121 typedef typename JacobianMatrix::block_type BlockType;
122 BlockType &J = jac[vIdx][vIdx];
123 for (int j = 0; j < BlockType::rows; ++j)
124 J[j][j] = 1.0;
125
126 // set residual for the vertex
127 res[vIdx] = 0;
128 }
129 }
130
131 auto applyDirichlet = [&] (const auto& scvI,
132 const auto& dirichletValues,
133 const auto eqIdx,
134 const auto pvIdx)
135 {
136 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
137
138 auto& row = jac[scvI.dofIndex()];
139 for (auto col = row.begin(); col != row.end(); ++col)
140 row[col.index()][eqIdx] = 0.0;
141
142 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
143
144 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
145 if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
146 {
147 const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
148 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
149 const auto end = jac[periodicDof].end();
150 for (auto it = jac[periodicDof].begin(); it != end; ++it)
151 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
152 }
153 };
154
155 this->asImp_().enforceDirichletConstraints(applyDirichlet);
156 }
157
162 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
163 {
164 this->asImp_().bindLocalViews();
165 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
166
167 auto applyDirichlet = [&] (const auto& scvI,
168 const auto& dirichletValues,
169 const auto eqIdx,
170 const auto pvIdx)
171 {
172 auto& row = jac[scvI.dofIndex()];
173 for (auto col = row.begin(); col != row.end(); ++col)
174 row[col.index()][eqIdx] = 0.0;
175
176 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
177 };
178
179 this->asImp_().enforceDirichletConstraints(applyDirichlet);
180 }
181
185 void assembleResidual(SolutionVector& res)
186 {
187 this->asImp_().bindLocalViews();
188 const auto residual = this->evalLocalResidual();
189
190 for (const auto& scv : scvs(this->fvGeometry()))
191 res[scv.dofIndex()] += residual[scv.localDofIndex()];
192
193 auto applyDirichlet = [&] (const auto& scvI,
194 const auto& dirichletValues,
195 const auto eqIdx,
196 const auto pvIdx)
197 {
198 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
199 };
200
201 this->asImp_().enforceDirichletConstraints(applyDirichlet);
202 }
203
205 template<typename ApplyFunction>
206 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
207 {
208 // enforce Dirichlet boundary conditions
209 this->asImp_().evalDirichletBoundaries(applyDirichlet);
210 // take care of internal Dirichlet constraints (if enabled)
211 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
212 }
213
217 template< typename ApplyDirichletFunctionType >
218 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
219 {
220 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
221 // and set the residual to (privar - dirichletvalue)
222 if (this->elemBcTypes().hasDirichlet())
223 {
224 for (const auto& scvI : scvs(this->fvGeometry()))
225 {
226 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
227 if (bcTypes.hasDirichlet())
228 {
229 const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
230
231 // set the Dirichlet conditions in residual and jacobian
232 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
233 {
234 if (bcTypes.isDirichlet(eqIdx))
235 {
236 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
237 assert(0 <= pvIdx && pvIdx < numEq);
238 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
239 }
240 }
241 }
242 }
243 }
244 }
245
246};
247
256template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true>
258
264template<class TypeTag, class Assembler>
265class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
266: public BoxLocalAssemblerBase<TypeTag, Assembler,
267 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true>
268{
276 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
277
280
282
283public:
284
286
293 template <class PartialReassembler = DefaultPartialReassembler>
294 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
295 const PartialReassembler* partialReassembler = nullptr)
296 {
297 // get some aliases for convenience
298 const auto& element = this->element();
299 const auto& fvGeometry = this->fvGeometry();
300 const auto& curSol = this->curSol();
301 auto&& curElemVolVars = this->curElemVolVars();
302
303 // get the vector of the actual element residuals
304 const auto origResiduals = this->evalLocalResidual();
305
307 // //
308 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
309 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
310 // actual element. In the actual element we evaluate the derivative of the entire residual. //
311 // //
313
314 // if all volvars in the stencil have to be updated or if it's enough to only update the
315 // volVars for the scv whose associated dof has been deflected
316 static const bool updateAllVolVars = getParamFromGroup<bool>(
317 this->problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
318 );
319
320 // create the element solution
321 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
322
323 // create the vector storing the partial derivatives
324 ElementResidualVector partialDerivs(fvGeometry.numScv());
325
326 Detail::VolVarsDeflectionHelper deflectionHelper(
327 [&] (const auto& scv) -> VolumeVariables& {
328 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
329 },
330 fvGeometry,
331 updateAllVolVars
332 );
333
334 // calculation of the derivatives
335 for (auto&& scv : scvs(fvGeometry))
336 {
337 // dof index and corresponding actual pri vars
338 const auto dofIdx = scv.dofIndex();
339 deflectionHelper.setCurrent(scv);
340
341 // calculate derivatives w.r.t to the privars at the dof at hand
342 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
343 {
344 partialDerivs = 0.0;
345
346 auto evalResiduals = [&](Scalar priVar)
347 {
348 // update the volume variables and compute element residual
349 elemSol[scv.localDofIndex()][pvIdx] = priVar;
350 deflectionHelper.deflect(elemSol, scv, this->problem());
351 return this->evalLocalResidual();
352 };
353
354 // derive the residuals numerically
355 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
356 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
357 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
358 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
359
360 // update the global stiffness matrix with the current partial derivatives
361 for (auto&& scvJ : scvs(fvGeometry))
362 {
363 // don't add derivatives for green vertices
364 if (!partialReassembler
365 || partialReassembler->vertexColor(scvJ.dofIndex()) != EntityColor::green)
366 {
367 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
368 {
369 // A[i][col][eqIdx][pvIdx] is the rate of change of
370 // the residual of equation 'eqIdx' at dof 'i'
371 // depending on the primary variable 'pvIdx' at dof
372 // 'col'.
373 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
374 }
375 }
376 }
377
378 // restore the original state of the scv's volume variables
379 deflectionHelper.restore(scv);
380
381 // restore the original element solution
382 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
383 }
384 }
385 return origResiduals;
386 }
387
388}; // implicit BoxAssembler with numeric Jacobian
389
395template<class TypeTag, class Assembler>
396class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
397: public BoxLocalAssemblerBase<TypeTag, Assembler,
398 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
399{
407 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
408
411
412public:
413
415
422 template <class PartialReassembler = DefaultPartialReassembler>
423 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
424 const PartialReassembler* partialReassembler = nullptr)
425 {
426 if (partialReassembler)
427 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
428
429 // get some aliases for convenience
430 const auto& element = this->element();
431 const auto& fvGeometry = this->fvGeometry();
432 const auto& curSol = this->curSol();
433 auto&& curElemVolVars = this->curElemVolVars();
434
435 // get the vecor of the actual element residuals
436 const auto origResiduals = this->evalLocalResidual();
437 const auto origStorageResiduals = this->evalLocalStorageResidual();
438
440 // //
441 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
442 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
443 // actual element. In the actual element we evaluate the derivative of the entire residual. //
444 // //
446
447 // create the element solution
448 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
449
450 // create the vector storing the partial derivatives
451 ElementResidualVector partialDerivs(fvGeometry.numScv());
452
453 // calculation of the derivatives
454 for (auto&& scv : scvs(fvGeometry))
455 {
456 // dof index and corresponding actual pri vars
457 const auto dofIdx = scv.dofIndex();
458 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
459 const VolumeVariables origVolVars(curVolVars);
460
461 // calculate derivatives w.r.t to the privars at the dof at hand
462 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
463 {
464 partialDerivs = 0.0;
465
466 auto evalStorage = [&](Scalar priVar)
467 {
468 // auto partialDerivsTmp = partialDerivs;
469 elemSol[scv.localDofIndex()][pvIdx] = priVar;
470 curVolVars.update(elemSol, this->problem(), element, scv);
471 return this->evalLocalStorageResidual();
472 };
473
474 // derive the residuals numerically
475 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
476 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
477 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
478 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
479
480 // update the global stiffness matrix with the current partial derivatives
481 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
482 {
483 // A[i][col][eqIdx][pvIdx] is the rate of change of
484 // the residual of equation 'eqIdx' at dof 'i'
485 // depending on the primary variable 'pvIdx' at dof
486 // 'col'.
487 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
488 }
489
490 // restore the original state of the scv's volume variables
491 curVolVars = origVolVars;
492
493 // restore the original element solution
494 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
495 // TODO additional dof dependencies
496 }
497 }
498 return origResiduals;
499 }
500}; // explicit BoxAssembler with numeric Jacobian
501
507template<class TypeTag, class Assembler>
508class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
509: public BoxLocalAssemblerBase<TypeTag, Assembler,
510 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
511{
517 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
518
519public:
520
522
529 template <class PartialReassembler = DefaultPartialReassembler>
530 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
531 const PartialReassembler* partialReassembler = nullptr)
532 {
533 if (partialReassembler)
534 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
535
536 // get some aliases for convenience
537 const auto& element = this->element();
538 const auto& fvGeometry = this->fvGeometry();
539 const auto& problem = this->problem();
540 const auto& curElemVolVars = this->curElemVolVars();
541 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
542
543 // get the vecor of the actual element residuals
544 const auto origResiduals = this->evalLocalResidual();
545
547 // //
548 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
549 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
550 // actual element. In the actual element we evaluate the derivative of the entire residual. //
551 // //
553
554 // calculation of the source and storage derivatives
555 for (const auto& scv : scvs(fvGeometry))
556 {
557 // dof index and corresponding actual pri vars
558 const auto dofIdx = scv.dofIndex();
559 const auto& volVars = curElemVolVars[scv];
560
561 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
562 // only if the problem is instationary we add derivative of storage term
563 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
564 if (!this->assembler().isStationaryProblem())
565 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
566 problem,
567 element,
568 fvGeometry,
569 volVars,
570 scv);
571
572 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
573 // add source term derivatives
574 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
575 problem,
576 element,
577 fvGeometry,
578 volVars,
579 scv);
580 }
581
582 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
583 for (const auto& scvf : scvfs(fvGeometry))
584 {
585 if (!scvf.boundary())
586 {
587 // add flux term derivatives
588 this->localResidual().addFluxDerivatives(A,
589 problem,
590 element,
591 fvGeometry,
592 curElemVolVars,
593 elemFluxVarsCache,
594 scvf);
595 }
596
597 // the boundary gets special treatment to simplify
598 // for the user
599 else
600 {
601 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
602 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
603 {
604 // add flux term derivatives
605 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
606 problem,
607 element,
608 fvGeometry,
609 curElemVolVars,
610 elemFluxVarsCache,
611 scvf);
612 }
613 }
614 }
615
616 return origResiduals;
617 }
618
619}; // implicit BoxAssembler with analytic Jacobian
620
626template<class TypeTag, class Assembler>
627class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
628: public BoxLocalAssemblerBase<TypeTag, Assembler,
629 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
630{
636 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
637
638public:
639
641
648 template <class PartialReassembler = DefaultPartialReassembler>
649 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
650 const PartialReassembler* partialReassembler = nullptr)
651 {
652 if (partialReassembler)
653 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
654
655 // get some aliases for convenience
656 const auto& element = this->element();
657 const auto& fvGeometry = this->fvGeometry();
658 const auto& problem = this->problem();
659 const auto& curElemVolVars = this->curElemVolVars();
660
661 // get the vecor of the actual element residuals
662 const auto origResiduals = this->evalLocalResidual();
663
665 // //
666 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
667 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
668 // actual element. In the actual element we evaluate the derivative of the entire residual. //
669 // //
671
672 // calculation of the source and storage derivatives
673 for (const auto& scv : scvs(fvGeometry))
674 {
675 // dof index and corresponding actual pri vars
676 const auto dofIdx = scv.dofIndex();
677 const auto& volVars = curElemVolVars[scv];
678
679 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
680 // only if the problem is instationary we add derivative of storage term
681 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
682 problem,
683 element,
684 fvGeometry,
685 volVars,
686 scv);
687 }
688
689 return origResiduals;
690 }
691
692}; // explicit BoxAssembler with analytic Jacobian
693
694} // end namespace Dumux
695
696#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
Detects which entries in the Jacobian have to be recomputed.
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
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==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
A base class for all local box assemblers.
Definition: boxlocalassembler.hh:57
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:162
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: boxlocalassembler.hh:218
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: boxlocalassembler.hh:206
void bindLocalViews()
Definition: boxlocalassembler.hh:70
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: boxlocalassembler.hh:185
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:82
An assembler for Jacobian and residual contribution per element (box methods)
Definition: boxlocalassembler.hh:257
Box local assembler using numeric differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:268
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:294
Box local assembler using numeric differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:399
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:423
Box local assembler using analytic differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:511
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:530
Box local assembler using analytic differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:630
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:649
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
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
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
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.
The local element solution class for the box method.