3.5-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
42namespace Dumux {
43
44#ifndef DOXYGEN
45
46namespace Detail {
47
48template<class VolVarAccessor, class FVElementGeometry>
49class VolVarsDeflectionHelper
50{
51 static constexpr int maxNumscvs = FVElementGeometry::maxNumElementScvs;
52
53 using SubControlVolume = typename FVElementGeometry::GridGeometry::SubControlVolume;
54 using VolumeVariablesRef = std::invoke_result_t<VolVarAccessor, const SubControlVolume&>;
55 using VolumeVariables = std::decay_t<VolumeVariablesRef>;
56 static_assert(std::is_lvalue_reference_v<VolumeVariablesRef>
57 && !std::is_const_v<std::remove_reference_t<VolumeVariablesRef>>);
58
59public:
60 VolVarsDeflectionHelper(VolVarAccessor&& accessor,
61 const FVElementGeometry& fvGeometry,
62 bool deflectAllVolVars)
63 : deflectAll_(deflectAllVolVars)
64 , accessor_(std::move(accessor))
65 , fvGeometry_(fvGeometry)
66 {
67 if (deflectAll_)
68 for (const auto& curScv : scvs(fvGeometry))
69 origVolVars_.push_back(accessor_(curScv));
70 }
71
72 void setCurrent(const SubControlVolume& scv)
73 {
74 if (!deflectAll_)
75 {
76 origVolVars_.clear();
77 origVolVars_.push_back(accessor_(scv));
78 }
79 }
80
81 template<class ElementSolution, class Problem>
82 void deflect(const ElementSolution& elemSol,
83 const SubControlVolume& scv,
84 const Problem& problem)
85 {
86 if (deflectAll_)
87 for (const auto& curScv : scvs(fvGeometry_))
88 accessor_(curScv).update(elemSol, problem, fvGeometry_.element(), curScv);
89 else
90 accessor_(scv).update(elemSol, problem, fvGeometry_.element(), scv);
91 }
92
93 void restore(const SubControlVolume& scv)
94 {
95 if (!deflectAll_)
96 accessor_(scv) = origVolVars_[0];
97 else
98 for (const auto& curScv : scvs(fvGeometry_))
99 accessor_(curScv) = origVolVars_[curScv.localDofIndex()];
100 }
101
102private:
103 const bool deflectAll_;
104 VolVarAccessor accessor_;
105 const FVElementGeometry& fvGeometry_;
106 Dune::ReservedVector<VolumeVariables, maxNumscvs> origVolVars_;
107};
108
109template<class Accessor, class FVElementGeometry>
110VolVarsDeflectionHelper(Accessor&&, FVElementGeometry&&) -> VolVarsDeflectionHelper<Accessor, std::decay_t<FVElementGeometry>>;
111
112} // end namespace Detail
113#endif // DOXYGEN
114
124template<class TypeTag, class Assembler, class Implementation, bool implicit>
125class BoxLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
126{
131 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
132
134
135public:
136
137 using ParentType::ParentType;
138
140 {
142 this->elemBcTypes().update(this->problem(), this->element(), this->fvGeometry());
143 }
144
145
150 template <class PartialReassembler = DefaultPartialReassembler>
151 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
152 const PartialReassembler* partialReassembler = nullptr)
153 {
154 this->asImp_().bindLocalViews();
155 const auto eIdxGlobal = this->assembler().gridGeometry().elementMapper().index(this->element());
156 if (partialReassembler
157 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
158 {
159 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
160 for (const auto& scv : scvs(this->fvGeometry()))
161 res[scv.dofIndex()] += residual[scv.localDofIndex()];
162 }
163 else if (!this->elementIsGhost())
164 {
165 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
166 for (const auto& scv : scvs(this->fvGeometry()))
167 res[scv.dofIndex()] += residual[scv.localDofIndex()];
168 }
169 else
170 {
171 using GridGeometry = typename GridVariables::GridGeometry;
172 using GridView = typename GridGeometry::GridView;
173 static constexpr auto dim = GridView::dimension;
174
175 int numVerticesLocal = this->element().subEntities(dim);
176
177 for (int i = 0; i < numVerticesLocal; ++i) {
178 auto vertex = this->element().template subEntity<dim>(i);
179
180 if (vertex.partitionType() == Dune::InteriorEntity ||
181 vertex.partitionType() == Dune::BorderEntity)
182 {
183 // do not change the non-ghost vertices
184 continue;
185 }
186
187 // set main diagonal entries for the vertex
188 int vIdx = this->assembler().gridGeometry().vertexMapper().index(vertex);
189
190 typedef typename JacobianMatrix::block_type BlockType;
191 BlockType &J = jac[vIdx][vIdx];
192 for (int j = 0; j < BlockType::rows; ++j)
193 J[j][j] = 1.0;
194
195 // set residual for the vertex
196 res[vIdx] = 0;
197 }
198 }
199
200 auto applyDirichlet = [&] (const auto& scvI,
201 const auto& dirichletValues,
202 const auto eqIdx,
203 const auto pvIdx)
204 {
205 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
206
207 auto& row = jac[scvI.dofIndex()];
208 for (auto col = row.begin(); col != row.end(); ++col)
209 row[col.index()][eqIdx] = 0.0;
210
211 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
212
213 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof
214 if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
215 {
216 const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
217 res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
218 const auto end = jac[periodicDof].end();
219 for (auto it = jac[periodicDof].begin(); it != end; ++it)
220 (*it) = periodicDof != it.index() ? 0.0 : 1.0;
221 }
222 };
223
224 this->asImp_().enforceDirichletConstraints(applyDirichlet);
225 }
226
231 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
232 {
233 this->asImp_().bindLocalViews();
234 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
235
236 auto applyDirichlet = [&] (const auto& scvI,
237 const auto& dirichletValues,
238 const auto eqIdx,
239 const auto pvIdx)
240 {
241 auto& row = jac[scvI.dofIndex()];
242 for (auto col = row.begin(); col != row.end(); ++col)
243 row[col.index()][eqIdx] = 0.0;
244
245 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
246 };
247
248 this->asImp_().enforceDirichletConstraints(applyDirichlet);
249 }
250
254 void assembleResidual(SolutionVector& res)
255 {
256 this->asImp_().bindLocalViews();
257 const auto residual = this->evalLocalResidual();
258
259 for (const auto& scv : scvs(this->fvGeometry()))
260 res[scv.dofIndex()] += residual[scv.localDofIndex()];
261
262 auto applyDirichlet = [&] (const auto& scvI,
263 const auto& dirichletValues,
264 const auto eqIdx,
265 const auto pvIdx)
266 {
267 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
268 };
269
270 this->asImp_().enforceDirichletConstraints(applyDirichlet);
271 }
272
274 template<typename ApplyFunction>
275 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
276 {
277 // enforce Dirichlet boundary conditions
278 this->asImp_().evalDirichletBoundaries(applyDirichlet);
279 // take care of internal Dirichlet constraints (if enabled)
280 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
281 }
282
286 template< typename ApplyDirichletFunctionType >
287 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
288 {
289 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
290 // and set the residual to (privar - dirichletvalue)
291 if (this->elemBcTypes().hasDirichlet())
292 {
293 for (const auto& scvI : scvs(this->fvGeometry()))
294 {
295 const auto bcTypes = this->elemBcTypes()[scvI.localDofIndex()];
296 if (bcTypes.hasDirichlet())
297 {
298 const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
299
300 // set the Dirichlet conditions in residual and jacobian
301 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
302 {
303 if (bcTypes.isDirichlet(eqIdx))
304 {
305 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
306 assert(0 <= pvIdx && pvIdx < numEq);
307 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
308 }
309 }
310 }
311 }
312 }
313 }
314
315};
316
325template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true>
327
333template<class TypeTag, class Assembler>
334class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
335: public BoxLocalAssemblerBase<TypeTag, Assembler,
336 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true>, true>
337{
345 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
346
349
350 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
351
352public:
353
355
362 template <class PartialReassembler = DefaultPartialReassembler>
363 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
364 const PartialReassembler* partialReassembler = nullptr)
365 {
366 // get some aliases for convenience
367 const auto& element = this->element();
368 const auto& fvGeometry = this->fvGeometry();
369 const auto& curSol = this->curSol();
370 auto&& curElemVolVars = this->curElemVolVars();
371
372 // get the vector of the actual element residuals
373 const auto origResiduals = this->evalLocalResidual();
374
376 // //
377 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
378 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
379 // actual element. In the actual element we evaluate the derivative of the entire residual. //
380 // //
382
383 // if all volvars in the stencil have to be updated or if it's enough to only update the
384 // volVars for the scv whose associated dof has been deflected
385 static const bool updateAllVolVars = getParamFromGroup<bool>(
386 this->problem().paramGroup(), "Assembly.BoxVolVarsDependOnAllElementDofs", false
387 );
388
389 // create the element solution
390 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
391
392 // create the vector storing the partial derivatives
393 ElementResidualVector partialDerivs(element.subEntities(dim));
394
395 Detail::VolVarsDeflectionHelper deflectionHelper(
396 [&] (const auto& scv) -> VolumeVariables& {
397 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
398 },
399 fvGeometry,
400 updateAllVolVars
401 );
402
403 // calculation of the derivatives
404 for (auto&& scv : scvs(fvGeometry))
405 {
406 // dof index and corresponding actual pri vars
407 const auto dofIdx = scv.dofIndex();
408 deflectionHelper.setCurrent(scv);
409
410 // calculate derivatives w.r.t to the privars at the dof at hand
411 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
412 {
413 partialDerivs = 0.0;
414
415 auto evalResiduals = [&](Scalar priVar)
416 {
417 // update the volume variables and compute element residual
418 elemSol[scv.localDofIndex()][pvIdx] = priVar;
419 deflectionHelper.deflect(elemSol, scv, this->problem());
420 return this->evalLocalResidual();
421 };
422
423 // derive the residuals numerically
424 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
425 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
426 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
427 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
428
429 // update the global stiffness matrix with the current partial derivatives
430 for (auto&& scvJ : scvs(fvGeometry))
431 {
432 // don't add derivatives for green vertices
433 if (!partialReassembler
434 || partialReassembler->vertexColor(scvJ.dofIndex()) != EntityColor::green)
435 {
436 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
437 {
438 // A[i][col][eqIdx][pvIdx] is the rate of change of
439 // the residual of equation 'eqIdx' at dof 'i'
440 // depending on the primary variable 'pvIdx' at dof
441 // 'col'.
442 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
443 }
444 }
445 }
446
447 // restore the original state of the scv's volume variables
448 deflectionHelper.restore(scv);
449
450 // restore the original element solution
451 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
452 }
453 }
454 return origResiduals;
455 }
456
457}; // implicit BoxAssembler with numeric Jacobian
458
464template<class TypeTag, class Assembler>
465class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
466: public BoxLocalAssemblerBase<TypeTag, Assembler,
467 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
468{
476 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
477
480
481public:
482
484
491 template <class PartialReassembler = DefaultPartialReassembler>
492 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
493 const PartialReassembler* partialReassembler = nullptr)
494 {
495 if (partialReassembler)
496 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
497
498 // get some aliases for convenience
499 const auto& element = this->element();
500 const auto& fvGeometry = this->fvGeometry();
501 const auto& curSol = this->curSol();
502 auto&& curElemVolVars = this->curElemVolVars();
503
504 // get the vecor of the acutal element residuals
505 const auto origResiduals = this->evalLocalResidual();
506 const auto origStorageResiduals = this->evalLocalStorageResidual();
507
509 // //
510 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
511 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
512 // actual element. In the actual element we evaluate the derivative of the entire residual. //
513 // //
515
516 // create the element solution
517 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
518
519 // create the vector storing the partial derivatives
520 ElementResidualVector partialDerivs(element.subEntities(dim));
521
522 // calculation of the derivatives
523 for (auto&& scv : scvs(fvGeometry))
524 {
525 // dof index and corresponding actual pri vars
526 const auto dofIdx = scv.dofIndex();
527 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
528 const VolumeVariables origVolVars(curVolVars);
529
530 // calculate derivatives w.r.t to the privars at the dof at hand
531 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
532 {
533 partialDerivs = 0.0;
534
535 auto evalStorage = [&](Scalar priVar)
536 {
537 // auto partialDerivsTmp = partialDerivs;
538 elemSol[scv.localDofIndex()][pvIdx] = priVar;
539 curVolVars.update(elemSol, this->problem(), element, scv);
540 return this->evalLocalStorageResidual();
541 };
542
543 // derive the residuals numerically
544 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
545 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
546 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
547 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
548
549 // update the global stiffness matrix with the current partial derivatives
550 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
551 {
552 // A[i][col][eqIdx][pvIdx] is the rate of change of
553 // the residual of equation 'eqIdx' at dof 'i'
554 // depending on the primary variable 'pvIdx' at dof
555 // 'col'.
556 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
557 }
558
559 // restore the original state of the scv's volume variables
560 curVolVars = origVolVars;
561
562 // restore the original element solution
563 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
564 // TODO additional dof dependencies
565 }
566 }
567 return origResiduals;
568 }
569}; // explicit BoxAssembler with numeric Jacobian
570
576template<class TypeTag, class Assembler>
577class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
578: public BoxLocalAssemblerBase<TypeTag, Assembler,
579 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
580{
586 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
587
588public:
589
591
598 template <class PartialReassembler = DefaultPartialReassembler>
599 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
600 const PartialReassembler* partialReassembler = nullptr)
601 {
602 if (partialReassembler)
603 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
604
605 // get some aliases for convenience
606 const auto& element = this->element();
607 const auto& fvGeometry = this->fvGeometry();
608 const auto& problem = this->problem();
609 const auto& curElemVolVars = this->curElemVolVars();
610 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
611
612 // get the vecor of the acutal element residuals
613 const auto origResiduals = this->evalLocalResidual();
614
616 // //
617 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
618 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
619 // actual element. In the actual element we evaluate the derivative of the entire residual. //
620 // //
622
623 // calculation of the source and storage derivatives
624 for (const auto& scv : scvs(fvGeometry))
625 {
626 // dof index and corresponding actual pri vars
627 const auto dofIdx = scv.dofIndex();
628 const auto& volVars = curElemVolVars[scv];
629
630 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
631 // only if the problem is instationary we add derivative of storage term
632 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
633 if (!this->assembler().isStationaryProblem())
634 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
635 problem,
636 element,
637 fvGeometry,
638 volVars,
639 scv);
640
641 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
642 // add source term derivatives
643 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
644 problem,
645 element,
646 fvGeometry,
647 volVars,
648 scv);
649 }
650
651 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
652 for (const auto& scvf : scvfs(fvGeometry))
653 {
654 if (!scvf.boundary())
655 {
656 // add flux term derivatives
657 this->localResidual().addFluxDerivatives(A,
658 problem,
659 element,
660 fvGeometry,
661 curElemVolVars,
662 elemFluxVarsCache,
663 scvf);
664 }
665
666 // the boundary gets special treatment to simplify
667 // for the user
668 else
669 {
670 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
671 if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
672 {
673 // add flux term derivatives
674 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
675 problem,
676 element,
677 fvGeometry,
678 curElemVolVars,
679 elemFluxVarsCache,
680 scvf);
681 }
682 }
683 }
684
685 return origResiduals;
686 }
687
688}; // implicit BoxAssembler with analytic Jacobian
689
695template<class TypeTag, class Assembler>
696class BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
697: public BoxLocalAssemblerBase<TypeTag, Assembler,
698 BoxLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
699{
705 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
706
707public:
708
710
717 template <class PartialReassembler = DefaultPartialReassembler>
718 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
719 const PartialReassembler* partialReassembler = nullptr)
720 {
721 if (partialReassembler)
722 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
723
724 // get some aliases for convenience
725 const auto& element = this->element();
726 const auto& fvGeometry = this->fvGeometry();
727 const auto& problem = this->problem();
728 const auto& curElemVolVars = this->curElemVolVars();
729
730 // get the vecor of the acutal element residuals
731 const auto origResiduals = this->evalLocalResidual();
732
734 // //
735 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
736 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
737 // actual element. In the actual element we evaluate the derivative of the entire residual. //
738 // //
740
741 // calculation of the source and storage derivatives
742 for (const auto& scv : scvs(fvGeometry))
743 {
744 // dof index and corresponding actual pri vars
745 const auto dofIdx = scv.dofIndex();
746 const auto& volVars = curElemVolVars[scv];
747
748 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
749 // only if the problem is instationary we add derivative of storage term
750 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
751 problem,
752 element,
753 fvGeometry,
754 volVars,
755 scv);
756 }
757
758 return origResiduals;
759 }
760
761}; // explicit BoxAssembler with analytic Jacobian
762
763} // end namespace Dumux
764
765#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==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
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
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:126
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:231
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: boxlocalassembler.hh:287
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: boxlocalassembler.hh:275
void bindLocalViews()
Definition: boxlocalassembler.hh:139
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: boxlocalassembler.hh:254
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:151
An assembler for Jacobian and residual contribution per element (box methods)
Definition: boxlocalassembler.hh:326
Box local assembler using numeric differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:337
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:363
Box local assembler using numeric differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:468
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:492
Box local assembler using analytic differentiation and implicit time discretization.
Definition: boxlocalassembler.hh:580
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:599
Box local assembler using analytic differentiation and explicit time discretization.
Definition: boxlocalassembler.hh:699
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:718
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 repect to a function parameter.
Definition: numericdifferentiation.hh:61
Declares all properties used in Dumux.
The local element solution class for the box method.