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