3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
subdomainboxlocalassembler.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 *****************************************************************************/
26#ifndef DUMUX_MULTIDOMAIN_BOX_LOCAL_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_BOX_LOCAL_ASSEMBLER_HH
28
29#include <dune/common/reservedvector.hh>
30#include <dune/common/indices.hh>
31#include <dune/common/hybridutilities.hh>
32#include <dune/grid/common/gridenums.hh> // for GhostEntity
33#include <dune/istl/matrixindexset.hh>
34
43
44namespace Dumux {
45
57template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool implicit>
58class SubDomainBoxLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
59{
61
63 using LocalResidualValues = GetPropType<TypeTag, Properties::NumEqVector>;
64 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
66 using SolutionVector = typename Assembler::SolutionVector;
69
71 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
72 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
73 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
74 using Scalar = typename GridVariables::Scalar;
75
76 using GridGeometry = typename GridVariables::GridGeometry;
77 using FVElementGeometry = typename GridGeometry::LocalView;
78 using SubControlVolume = typename GridGeometry::SubControlVolume;
79 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
80 using Extrusion = Extrusion_t<GridGeometry>;
81 using GridView = typename GridGeometry::GridView;
82 using Element = typename GridView::template Codim<0>::Entity;
83
84 using CouplingManager = typename Assembler::CouplingManager;
85
86 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
87
88public:
90 static constexpr auto domainId = typename Dune::index_constant<id>();
92 using ParentType::ParentType;
93
94 // the constructor
95 explicit SubDomainBoxLocalAssemblerBase(const Assembler& assembler,
96 const Element& element,
97 const SolutionVector& curSol,
98 CouplingManager& couplingManager)
100 element,
101 curSol,
102 localView(assembler.gridGeometry(domainId)),
103 localView(assembler.gridVariables(domainId).curGridVolVars()),
104 localView(assembler.gridVariables(domainId).prevGridVolVars()),
105 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
107 (element.partitionType() == Dune::GhostEntity))
108 , couplingManager_(couplingManager)
109 {}
110
115 template<class JacobianMatrixRow, class GridVariablesTuple>
116 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSolutionVector& res, GridVariablesTuple& gridVariables)
117 {
118 this->asImp_().bindLocalViews();
119 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
120
121 if (!this->elementIsGhost())
122 {
123 // for the diagonal jacobian block
124 // forward to the internal implementation
125 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jacRow[domainId], *std::get<domainId>(gridVariables));
126
127 // update the residual vector
128 for (const auto& scv : scvs(this->fvGeometry()))
129 res[scv.dofIndex()] += residual[scv.localDofIndex()];
130
131 // assemble the coupling blocks
132 using namespace Dune::Hybrid;
133 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
134 {
135 if (i != id)
136 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
137 });
138 }
139 else
140 {
141 using GridGeometry = typename GridVariables::GridGeometry;
142 using GridView = typename GridGeometry::GridView;
143 static constexpr auto dim = GridView::dimension;
144
145 int numVerticesLocal = this->element().subEntities(dim);
146
147 for (int i = 0; i < numVerticesLocal; ++i)
148 {
149 const auto vertex = this->element().template subEntity<dim>(i);
150
151 if (vertex.partitionType() == Dune::InteriorEntity ||
152 vertex.partitionType() == Dune::BorderEntity)
153 {
154 // do not change the non-ghost vertices
155 continue;
156 }
157
158 // set main diagonal entries for the vertex
159 const auto vIdx = this->assembler().gridGeometry(domainId).vertexMapper().index(vertex);
160
161 typedef typename JacobianMatrix::block_type BlockType;
162 BlockType &J = jacRow[domainId][vIdx][vIdx];
163 for (int j = 0; j < BlockType::rows; ++j)
164 J[j][j] = 1.0;
165
166 // set residual for the vertex
167 res[vIdx] = 0;
168 }
169 }
170
171 // lambda for the incorporation of Dirichlet Bcs
172 auto applyDirichlet = [&] (const auto& scvI,
173 const auto& dirichletValues,
174 const auto eqIdx,
175 const auto pvIdx)
176 {
177 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
178
179 // in explicit schemes we only have entries on the diagonal
180 // and thus don't have to do anything with off-diagonal entries
181 if (implicit)
182 {
183 for (const auto& scvJ : scvs(this->fvGeometry()))
184 jacRow[domainId][scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
185 }
186
187 jacRow[domainId][scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
188 };
189
190 // incorporate Dirichlet BCs
191 this->asImp_().enforceDirichletConstraints(applyDirichlet);
192 }
193
198 template<std::size_t otherId, class JacRow, class GridVariables,
199 typename std::enable_if_t<(otherId == id), int> = 0>
200 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
201 const ElementResidualVector& res, GridVariables& gridVariables)
202 {}
203
207 template<std::size_t otherId, class JacRow, class GridVariables,
208 typename std::enable_if_t<(otherId != id), int> = 0>
209 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
210 const ElementResidualVector& res, GridVariables& gridVariables)
211 {
212 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
213 }
214
218 void assembleResidual(SubSolutionVector& res)
219 {
220 this->asImp_().bindLocalViews();
221 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
222
223 const auto residual = this->evalLocalResidual();
224 for (const auto& scv : scvs(this->fvGeometry()))
225 res[scv.dofIndex()] += residual[scv.localDofIndex()];
226
227 auto applyDirichlet = [&] (const auto& scvI,
228 const auto& dirichletValues,
229 const auto eqIdx,
230 const auto pvIdx)
231 {
232 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
233 };
234
235 this->asImp_().enforceDirichletConstraints(applyDirichlet);
236 }
237
241 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
242 {
243 // initialize the residual vector for all scvs in this element
244 ElementResidualVector residual(this->fvGeometry().numScv());
245
246 // evaluate the volume terms (storage + source terms)
247 // forward to the local residual specialized for the discretization methods
248 for (auto&& scv : scvs(this->fvGeometry()))
249 {
250 const auto& curVolVars = elemVolVars[scv];
251 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
252 source *= -Extrusion::volume(scv)*curVolVars.extrusionFactor();
253 residual[scv.localDofIndex()] = std::move(source);
254 }
255
256 return residual;
257 }
258
262 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
263 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
264
266 template<typename ApplyFunction>
267 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
268 {
269 // enforce Dirichlet boundary conditions
270 this->asImp_().evalDirichletBoundaries(applyDirichlet);
271 // take care of internal Dirichlet constraints (if enabled)
272 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
273 }
274
279 template< typename ApplyDirichletFunctionType >
280 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
281 {
282 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
283 // and set the residual to (privar - dirichletvalue)
284 if (this->elemBcTypes().hasDirichlet())
285 {
286 for (const auto& scvI : scvs(this->fvGeometry()))
287 {
288 const auto bcTypes = this->elemBcTypes()[scvI.localDofIndex()];
289 if (bcTypes.hasDirichlet())
290 {
291 const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
292
293 // set the dirichlet conditions in residual and jacobian
294 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
295 {
296 if (bcTypes.isDirichlet(eqIdx))
297 {
298 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
299 assert(0 <= pvIdx && pvIdx < numEq);
300 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
301 }
302 }
303 }
304 }
305 }
306 }
307
312 {
313 // get some references for convenience
314 const auto& element = this->element();
315 const auto& curSol = this->curSol()[domainId];
316 auto&& fvGeometry = this->fvGeometry();
317 auto&& curElemVolVars = this->curElemVolVars();
318 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
319
320 // bind the caches
321 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
322 fvGeometry.bind(element);
323
324 if (implicit)
325 {
328 if (!this->assembler().isStationaryProblem())
329 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
330 }
331 else
332 {
333 auto& prevElemVolVars = this->prevElemVolVars();
334 const auto& prevSol = this->assembler().prevSol()[domainId];
335
337 prevElemVolVars.bind(element, fvGeometry, prevSol);
339 }
340 }
341
343 const Problem& problem() const
344 { return this->assembler().problem(domainId); }
345
347 CouplingManager& couplingManager()
348 { return couplingManager_; }
349
350private:
351 CouplingManager& couplingManager_;
352};
353
365template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
367
374template<std::size_t id, class TypeTag, class Assembler>
375class SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
376: public SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler,
378{
379 using ThisType = SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
380 using ParentType = SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
383 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
384
385 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
387 using FVElementGeometry = typename GridGeometry::LocalView;
388 using GridView = typename GridGeometry::GridView;
389 using Element = typename GridView::template Codim<0>::Entity;
390
392 enum { dim = GridView::dimension };
393
394 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
395 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
396 static constexpr auto domainI = Dune::index_constant<id>();
397
398public:
400
407 template<class JacobianMatrixDiagBlock, class GridVariables>
408 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
409 {
411 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
412 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
413 // actual element. In the actual element we evaluate the derivative of the entire residual. //
415
416 // get the vecor of the acutal element residuals
417 const auto origResiduals = this->evalLocalResidual();
418
420 // compute the derivatives of this element with respect to all of the element's dofs and add them to the Jacobian //
422 const auto& element = this->element();
423 const auto& fvGeometry = this->fvGeometry();
424 const auto& curSol = this->curSol()[domainI];
425 auto&& curElemVolVars = this->curElemVolVars();
426
427 // create the element solution
428 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
429
430 auto partialDerivs = origResiduals;
431 partialDerivs = 0.0;
432
433 for (auto&& scv : scvs(fvGeometry))
434 {
435 // dof index and corresponding actual pri vars
436 const auto dofIdx = scv.dofIndex();
437 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
438 const VolumeVariables origVolVars(curVolVars);
439
440 // calculate derivatives w.r.t to the privars at the dof at hand
441 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
442 {
443 partialDerivs = 0.0;
444
445 auto evalResiduals = [&](Scalar priVar)
446 {
447 // update the volume variables and compute element residual
448 const auto localDofIndex = scv.localDofIndex();
449 elemSol[localDofIndex][pvIdx] = priVar;
450 curVolVars.update(elemSol, this->problem(), element, scv);
451 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[localDofIndex], pvIdx);
452 return this->evalLocalResidual();
453 };
454
455 // derive the residuals numerically
456 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
457 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
458 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
459 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
460
461 // update the global stiffness matrix with the current partial derivatives
462 for (auto&& scvJ : scvs(fvGeometry))
463 {
464 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
465 {
466 // A[i][col][eqIdx][pvIdx] is the rate of change of
467 // the residual of equation 'eqIdx' at dof 'i'
468 // depending on the primary variable 'pvIdx' at dof
469 // 'col'.
470 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
471 }
472 }
473
474 // restore the original state of the scv's volume variables
475 curVolVars = origVolVars;
476
477 // restore the original element solution and coupling context
478 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
479 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
480 }
481 }
482
483 // evaluate additional derivatives that might arise from the coupling
484 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
485
486 return origResiduals;
487 }
488
495 template<std::size_t otherId, class JacobianBlock, class GridVariables>
496 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
497 const ElementResidualVector& res, GridVariables& gridVariables)
498 {
500 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
502
503 // get some aliases for convenience
504 const auto& element = this->element();
505 const auto& fvGeometry = this->fvGeometry();
506 auto&& curElemVolVars = this->curElemVolVars();
507 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
508
509 // get element stencil informations
510 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
511
512 // convenience lambda for call to update self
513 auto updateCoupledVariables = [&] ()
514 {
515 // Update ourself after the context has been modified. Depending on the
516 // type of caching, other objects might have to be updated. All ifs can be optimized away.
517 if (enableGridFluxVarsCache)
518 {
519 if (enableGridVolVarsCache)
520 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
521 else
522 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
523 }
524 else
525 {
526 if (enableGridVolVarsCache)
527 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
528 else
529 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
530 }
531 };
532
533 const auto& curSolJ = this->curSol()[domainJ];
534 for (const auto globalJ : stencil)
535 {
536 // undeflected privars and privars to be deflected
537 const auto origPriVarsJ = curSolJ[globalJ];
538 auto priVarsJ = origPriVarsJ;
539
540 // the undeflected coupling residual
541 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
542
543 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
544 {
545 auto evalCouplingResidual = [&](Scalar priVar)
546 {
547 priVarsJ[pvIdx] = priVar;
548 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
549 updateCoupledVariables();
550 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
551 };
552
553 // derive the residuals numerically
554 ElementResidualVector partialDerivs(element.subEntities(dim));
555
556 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
557 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
558 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
559
560 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
561 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
562
563 // update the global stiffness matrix with the current partial derivatives
564 for (auto&& scv : scvs(fvGeometry))
565 {
566 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
567 {
568 // A[i][col][eqIdx][pvIdx] is the rate of change of
569 // the residual of equation 'eqIdx' at dof 'i'
570 // depending on the primary variable 'pvIdx' at dof
571 // 'col'.
572
573 // If the dof is coupled by a Dirichlet condition,
574 // set the derived value only once (i.e. overwrite existing values).
575 // For other dofs, add the contribution of the partial derivative.
576 const auto bcTypes = this->elemBcTypes()[scv.localDofIndex()];
577 if (bcTypes.isCouplingDirichlet(eqIdx))
578 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
579 else if (bcTypes.isDirichlet(eqIdx))
580 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
581 else
582 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
583 }
584 }
585
586 // restore the current element solution
587 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
588
589 // restore the undeflected state of the coupling context
590 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
591 }
592
593 // Restore original state of the flux vars cache and/or vol vars.
594 // This has to be done in case they depend on variables of domainJ before
595 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
596 // the next "origResidual" will be incorrect.
597 updateCoupledVariables();
598 }
599 }
600};
601
608template<std::size_t id, class TypeTag, class Assembler>
609class SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
610: public SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler,
611 SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false >
612{
613 using ThisType = SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
614 using ParentType = SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/false>;
617 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
618 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
620 using FVElementGeometry = typename GridGeometry::LocalView;
621 using GridView = typename GridGeometry::GridView;
622 using Element = typename GridView::template Codim<0>::Entity;
623
625 enum { dim = GridView::dimension };
626
627 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
628 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
629 static constexpr auto domainI = Dune::index_constant<id>();
630
631public:
633
640 template<class JacobianMatrixDiagBlock, class GridVariables>
641 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
642 {
643 if (this->assembler().isStationaryProblem())
644 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
645
646 // get some aliases for convenience
647 const auto& element = this->element();
648 const auto& fvGeometry = this->fvGeometry();
649 const auto& curSol = this->curSol()[domainI];
650 auto&& curElemVolVars = this->curElemVolVars();
651
652 // get the vector of the actual (undeflected) element residuals
653 const auto origResiduals = this->evalLocalResidual();
654 const auto origStorageResiduals = this->evalLocalStorageResidual();
655
657 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
658 // neighboring elements all derivatives are zero. For the assembled element only the storage //
659 // derivatives are non-zero. //
661
662 // create the element solution
663 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
664
665 // create the vector storing the partial derivatives
666 ElementResidualVector partialDerivs(element.subEntities(dim));
667
668 // calculation of the derivatives
669 for (auto&& scv : scvs(fvGeometry))
670 {
671 // dof index and corresponding actual pri vars
672 const auto dofIdx = scv.dofIndex();
673 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
674 const VolumeVariables origVolVars(curVolVars);
675
676 // calculate derivatives w.r.t to the privars at the dof at hand
677 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
678 {
679 partialDerivs = 0.0;
680
681 auto evalStorage = [&](Scalar priVar)
682 {
683 // auto partialDerivsTmp = partialDerivs;
684 elemSol[scv.localDofIndex()][pvIdx] = priVar;
685 curVolVars.update(elemSol, this->problem(), element, scv);
686 return this->evalLocalStorageResidual();
687 };
688
689 // derive the residuals numerically
690 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
691 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
692 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
693 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
694
695 // update the global stiffness matrix with the current partial derivatives
696 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
697 {
698 // A[i][col][eqIdx][pvIdx] is the rate of change of
699 // the residual of equation 'eqIdx' at dof 'i'
700 // depending on the primary variable 'pvIdx' at dof
701 // 'col'.
702 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
703 }
704
705 // restore the original state of the scv's volume variables
706 curVolVars = origVolVars;
707
708 // restore the original element solution
709 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
710 // TODO additional dof dependencies
711 }
712 }
713
714 // return the undeflected residual
715 return origResiduals;
716 }
717
724 template<std::size_t otherId, class JacobianBlock, class GridVariables>
725 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
726 const ElementResidualVector& res, GridVariables& gridVariables)
727 {}
728};
729
730} // end namespace Dumux
731
732#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.
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.
Helper classes to compute the integration elements.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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
Definition: common/pdesolver.hh:35
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:269
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:245
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:273
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:277
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:257
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
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
Definition: multidomain/couplingmanager.hh:46
A base class for all box local assemblers.
Definition: subdomainboxlocalassembler.hh:59
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: subdomainboxlocalassembler.hh:267
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainboxlocalassembler.hh:311
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainboxlocalassembler.hh:90
ElementResidualVector evalLocalSourceResidual(const Element &element, const ElementVolumeVariables &elemVolVars) const
Evaluates the local source term for an element and given element volume variables.
Definition: subdomainboxlocalassembler.hh:241
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainboxlocalassembler.hh:347
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainboxlocalassembler.hh:262
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Incorporate Dirichlet boundary conditions.
Definition: subdomainboxlocalassembler.hh:280
void assembleResidual(SubSolutionVector &res)
Assemble the residual vector entries only.
Definition: subdomainboxlocalassembler.hh:218
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacRow &jacRow, const ElementResidualVector &res, GridVariables &gridVariables)
Assemble the entries in a coupling block of the jacobian. There is no coupling block between a domain...
Definition: subdomainboxlocalassembler.hh:200
const Problem & problem() const
return reference to the underlying problem
Definition: subdomainboxlocalassembler.hh:343
void assembleJacobianAndResidual(JacobianMatrixRow &jacRow, SubSolutionVector &res, GridVariablesTuple &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: subdomainboxlocalassembler.hh:116
SubDomainBoxLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainboxlocalassembler.hh:95
The box scheme multidomain local assembler.
Definition: subdomainboxlocalassembler.hh:366
Box scheme multi domain local assembler using numeric differentiation and implicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:378
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:408
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:496
Box scheme multi domain local assembler using numeric differentiation and explicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:612
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: subdomainboxlocalassembler.hh:641
void assembleJacobianCoupling(Dune::index_constant< otherId > domainJ, JacobianBlock &A, const ElementResidualVector &res, GridVariables &gridVariables)
Computes the coupling derivatives with respect to the given element and adds them to the global matri...
Definition: subdomainboxlocalassembler.hh:725
Declares all properties used in Dumux.