3.2-git
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
42
43namespace Dumux {
44
56template<std::size_t id, class TypeTag, class Assembler, class Implementation, bool implicit>
57class SubDomainBoxLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
58{
60
62 using LocalResidualValues = GetPropType<TypeTag, Properties::NumEqVector>;
63 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
65 using SolutionVector = typename Assembler::SolutionVector;
68
70 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
71 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
72 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
73 using Scalar = typename GridVariables::Scalar;
74
75 using GridGeometry = typename GridVariables::GridGeometry;
76 using FVElementGeometry = typename GridGeometry::LocalView;
77 using SubControlVolume = typename GridGeometry::SubControlVolume;
78 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
79 using GridView = typename GridGeometry::GridView;
80 using Element = typename GridView::template Codim<0>::Entity;
81
82 using CouplingManager = typename Assembler::CouplingManager;
83
84 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
85
86public:
88 static constexpr auto domainId = typename Dune::index_constant<id>();
90 using ParentType::ParentType;
91
92 // the constructor
93 explicit SubDomainBoxLocalAssemblerBase(const Assembler& assembler,
94 const Element& element,
95 const SolutionVector& curSol,
96 CouplingManager& couplingManager)
98 element,
99 curSol,
100 localView(assembler.gridGeometry(domainId)),
101 localView(assembler.gridVariables(domainId).curGridVolVars()),
102 localView(assembler.gridVariables(domainId).prevGridVolVars()),
103 localView(assembler.gridVariables(domainId).gridFluxVarsCache()),
105 (element.partitionType() == Dune::GhostEntity))
106 , couplingManager_(couplingManager)
107 {}
108
113 template<class JacobianMatrixRow, class GridVariablesTuple>
114 void assembleJacobianAndResidual(JacobianMatrixRow& jacRow, SubSolutionVector& res, GridVariablesTuple& gridVariables)
115 {
116 this->asImp_().bindLocalViews();
117 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
118
119 // for the diagonal jacobian block
120 // forward to the internal implementation
121 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jacRow[domainId], *std::get<domainId>(gridVariables));
122
123 // update the residual vector
124 for (const auto& scv : scvs(this->fvGeometry()))
125 res[scv.dofIndex()] += residual[scv.localDofIndex()];
126
127 // assemble the coupling blocks
128 using namespace Dune::Hybrid;
129 forEach(integralRange(Dune::Hybrid::size(jacRow)), [&](auto&& i)
130 {
131 if (i != id)
132 this->assembleJacobianCoupling(i, jacRow, residual, gridVariables);
133 });
134
135 // lambda for the incorporation of Dirichlet Bcs
136 auto applyDirichlet = [&] (const auto& scvI,
137 const auto& dirichletValues,
138 const auto eqIdx,
139 const auto pvIdx)
140 {
141 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
142
143 // in explicit schemes we only have entries on the diagonal
144 // and thus don't have to do anything with off-diagonal entries
145 if (implicit)
146 {
147 for (const auto& scvJ : scvs(this->fvGeometry()))
148 jacRow[domainId][scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
149 }
150
151 jacRow[domainId][scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
152 };
153
154 // incorporate Dirichlet BCs
155 this->asImp_().enforceDirichletConstraints(applyDirichlet);
156 }
157
162 template<std::size_t otherId, class JacRow, class GridVariables,
163 typename std::enable_if_t<(otherId == id), int> = 0>
164 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
165 const ElementResidualVector& res, GridVariables& gridVariables)
166 {}
167
171 template<std::size_t otherId, class JacRow, class GridVariables,
172 typename std::enable_if_t<(otherId != id), int> = 0>
173 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacRow& jacRow,
174 const ElementResidualVector& res, GridVariables& gridVariables)
175 {
176 this->asImp_().assembleJacobianCoupling(domainJ, jacRow[domainJ], res, *std::get<domainId>(gridVariables));
177 }
178
182 void assembleResidual(SubSolutionVector& res)
183 {
184 this->asImp_().bindLocalViews();
185 this->elemBcTypes().update(problem(), this->element(), this->fvGeometry());
186
187 const auto residual = this->evalLocalResidual();
188 for (const auto& scv : scvs(this->fvGeometry()))
189 res[scv.dofIndex()] += residual[scv.localDofIndex()];
190
191 auto applyDirichlet = [&] (const auto& scvI,
192 const auto& dirichletValues,
193 const auto eqIdx,
194 const auto pvIdx)
195 {
196 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
197 };
198
199 this->asImp_().enforceDirichletConstraints(applyDirichlet);
200 }
201
205 ElementResidualVector evalLocalSourceResidual(const Element& element, const ElementVolumeVariables& elemVolVars) const
206 {
207 // initialize the residual vector for all scvs in this element
208 ElementResidualVector residual(this->fvGeometry().numScv());
209
210 // evaluate the volume terms (storage + source terms)
211 // forward to the local residual specialized for the discretization methods
212 for (auto&& scv : scvs(this->fvGeometry()))
213 {
214 const auto& curVolVars = elemVolVars[scv];
215 auto source = this->localResidual().computeSource(problem(), element, this->fvGeometry(), elemVolVars, scv);
216 source *= -scv.volume()*curVolVars.extrusionFactor();
217 residual[scv.localDofIndex()] = std::move(source);
218 }
219
220 return residual;
221 }
222
226 ElementResidualVector evalLocalSourceResidual(const Element& neighbor) const
227 { return this->evalLocalSourceResidual(neighbor, implicit ? this->curElemVolVars() : this->prevElemVolVars()); }
228
230 template<typename ApplyFunction>
231 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
232 {
233 // enforce Dirichlet boundary conditions
234 this->asImp_().evalDirichletBoundaries(applyDirichlet);
235 // take care of internal Dirichlet constraints (if enabled)
236 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
237 }
238
243 template< typename ApplyDirichletFunctionType >
244 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
245 {
246 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
247 // and set the residual to (privar - dirichletvalue)
248 if (this->elemBcTypes().hasDirichlet())
249 {
250 for (const auto& scvI : scvs(this->fvGeometry()))
251 {
252 const auto bcTypes = this->elemBcTypes()[scvI.localDofIndex()];
253 if (bcTypes.hasDirichlet())
254 {
255 const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
256
257 // set the dirichlet conditions in residual and jacobian
258 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
259 {
260 if (bcTypes.isDirichlet(eqIdx))
261 {
262 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
263 assert(0 <= pvIdx && pvIdx < numEq);
264 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
265 }
266 }
267 }
268 }
269 }
270 }
271
276 {
277 // get some references for convenience
278 const auto& element = this->element();
279 const auto& curSol = this->curSol()[domainId];
280 auto&& fvGeometry = this->fvGeometry();
281 auto&& curElemVolVars = this->curElemVolVars();
282 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
283
284 // bind the caches
285 couplingManager_.bindCouplingContext(domainId, element, this->assembler());
286 fvGeometry.bind(element);
287
288 if (implicit)
289 {
292 if (!this->assembler().isStationaryProblem())
293 this->prevElemVolVars().bindElement(element, fvGeometry, this->assembler().prevSol()[domainId]);
294 }
295 else
296 {
297 auto& prevElemVolVars = this->prevElemVolVars();
298 const auto& prevSol = this->assembler().prevSol()[domainId];
299
301 prevElemVolVars.bind(element, fvGeometry, prevSol);
303 }
304 }
305
307 const Problem& problem() const
308 { return this->assembler().problem(domainId); }
309
311 CouplingManager& couplingManager()
312 { return couplingManager_; }
313
314private:
315 CouplingManager& couplingManager_;
316};
317
329template<std::size_t id, class TypeTag, class Assembler, DiffMethod DM = DiffMethod::numeric, bool implicit = true>
331
338template<std::size_t id, class TypeTag, class Assembler>
339class SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>
340: public SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler,
342{
343 using ThisType = SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true>;
344 using ParentType = SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/true>;
347 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
348
349 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
351 using FVElementGeometry = typename GridGeometry::LocalView;
352 using GridView = typename GridGeometry::GridView;
353 using Element = typename GridView::template Codim<0>::Entity;
354
356 enum { dim = GridView::dimension };
357
358 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
359 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
360 static constexpr auto domainI = Dune::index_constant<id>();
361
362public:
364
371 template<class JacobianMatrixDiagBlock, class GridVariables>
372 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
373 {
375 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
376 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
377 // actual element. In the actual element we evaluate the derivative of the entire residual. //
379
380 // get the vecor of the acutal element residuals
381 const auto origResiduals = this->evalLocalResidual();
382
384 // compute the derivatives of this element with respect to all of the element's dofs and add them to the Jacobian //
386 const auto& element = this->element();
387 const auto& fvGeometry = this->fvGeometry();
388 const auto& curSol = this->curSol()[domainI];
389 auto&& curElemVolVars = this->curElemVolVars();
390
391 // create the element solution
392 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
393
394 auto partialDerivs = origResiduals;
395 partialDerivs = 0.0;
396
397 for (auto&& scv : scvs(fvGeometry))
398 {
399 // dof index and corresponding actual pri vars
400 const auto dofIdx = scv.dofIndex();
401 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
402 const VolumeVariables origVolVars(curVolVars);
403
404 // calculate derivatives w.r.t to the privars at the dof at hand
405 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
406 {
407 partialDerivs = 0.0;
408
409 auto evalResiduals = [&](Scalar priVar)
410 {
411 // update the volume variables and compute element residual
412 const auto localDofIndex = scv.localDofIndex();
413 elemSol[localDofIndex][pvIdx] = priVar;
414 curVolVars.update(elemSol, this->problem(), element, scv);
415 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[localDofIndex], pvIdx);
416 return this->evalLocalResidual();
417 };
418
419 // derive the residuals numerically
420 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
421 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
422 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
423 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
424
425 // update the global stiffness matrix with the current partial derivatives
426 for (auto&& scvJ : scvs(fvGeometry))
427 {
428 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
429 {
430 // A[i][col][eqIdx][pvIdx] is the rate of change of
431 // the residual of equation 'eqIdx' at dof 'i'
432 // depending on the primary variable 'pvIdx' at dof
433 // 'col'.
434 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
435 }
436 }
437
438 // restore the original state of the scv's volume variables
439 curVolVars = origVolVars;
440
441 // restore the original element solution and coupling context
442 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
443 this->couplingManager().updateCouplingContext(domainI, *this, domainI, scv.dofIndex(), elemSol[scv.localDofIndex()], pvIdx);
444 }
445 }
446
447 // evaluate additional derivatives that might arise from the coupling
448 this->couplingManager().evalAdditionalDomainDerivatives(domainI, *this, origResiduals, A, gridVariables);
449
450 return origResiduals;
451 }
452
459 template<std::size_t otherId, class JacobianBlock, class GridVariables>
460 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
461 const ElementResidualVector& res, GridVariables& gridVariables)
462 {
464 // Calculate derivatives of all dofs in the element with respect to all dofs in the coupling stencil. //
466
467 // get some aliases for convenience
468 const auto& element = this->element();
469 const auto& fvGeometry = this->fvGeometry();
470 auto&& curElemVolVars = this->curElemVolVars();
471 auto&& elemFluxVarsCache = this->elemFluxVarsCache();
472
473 // get element stencil informations
474 const auto& stencil = this->couplingManager().couplingStencil(domainI, element, domainJ);
475
476 // convenience lambda for call to update self
477 auto updateCoupledVariables = [&] ()
478 {
479 // Update ourself after the context has been modified. Depending on the
480 // type of caching, other objects might have to be updated. All ifs can be optimized away.
481 if (enableGridFluxVarsCache)
482 {
483 if (enableGridVolVarsCache)
484 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), gridVariables.gridFluxVarsCache());
485 else
486 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, gridVariables.gridFluxVarsCache());
487 }
488 else
489 {
490 if (enableGridVolVarsCache)
491 this->couplingManager().updateCoupledVariables(domainI, *this, gridVariables.curGridVolVars(), elemFluxVarsCache);
492 else
493 this->couplingManager().updateCoupledVariables(domainI, *this, curElemVolVars, elemFluxVarsCache);
494 }
495 };
496
497 const auto& curSolJ = this->curSol()[domainJ];
498 for (const auto globalJ : stencil)
499 {
500 // undeflected privars and privars to be deflected
501 const auto origPriVarsJ = curSolJ[globalJ];
502 auto priVarsJ = origPriVarsJ;
503
504 // the undeflected coupling residual
505 const auto origResidual = this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
506
507 for (int pvIdx = 0; pvIdx < JacobianBlock::block_type::cols; ++pvIdx)
508 {
509 auto evalCouplingResidual = [&](Scalar priVar)
510 {
511 priVarsJ[pvIdx] = priVar;
512 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
513 updateCoupledVariables();
514 return this->couplingManager().evalCouplingResidual(domainI, *this, domainJ, globalJ);
515 };
516
517 // derive the residuals numerically
518 ElementResidualVector partialDerivs(element.subEntities(dim));
519
520 const auto& paramGroup = this->assembler().problem(domainJ).paramGroup();
521 static const int numDiffMethod = getParamFromGroup<int>(paramGroup, "Assembly.NumericDifferenceMethod");
522 static const auto epsCoupl = this->couplingManager().numericEpsilon(domainJ, paramGroup);
523
524 NumericDifferentiation::partialDerivative(evalCouplingResidual, origPriVarsJ[pvIdx], partialDerivs, origResidual,
525 epsCoupl(origPriVarsJ[pvIdx], pvIdx), numDiffMethod);
526
527 // update the global stiffness matrix with the current partial derivatives
528 for (auto&& scv : scvs(fvGeometry))
529 {
530 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
531 {
532 // A[i][col][eqIdx][pvIdx] is the rate of change of
533 // the residual of equation 'eqIdx' at dof 'i'
534 // depending on the primary variable 'pvIdx' at dof
535 // 'col'.
536
537 // If the dof is coupled by a Dirichlet condition,
538 // set the derived value only once (i.e. overwrite existing values).
539 // For other dofs, add the contribution of the partial derivative.
540 const auto bcTypes = this->elemBcTypes()[scv.localDofIndex()];
541 if (bcTypes.isCouplingDirichlet(eqIdx))
542 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = partialDerivs[scv.localDofIndex()][eqIdx];
543 else if (bcTypes.isDirichlet(eqIdx))
544 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] = 0.0;
545 else
546 A[scv.dofIndex()][globalJ][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
547 }
548 }
549
550 // restore the current element solution
551 priVarsJ[pvIdx] = origPriVarsJ[pvIdx];
552
553 // restore the undeflected state of the coupling context
554 this->couplingManager().updateCouplingContext(domainI, *this, domainJ, globalJ, priVarsJ, pvIdx);
555 }
556
557 // Restore original state of the flux vars cache and/or vol vars.
558 // This has to be done in case they depend on variables of domainJ before
559 // we continue with the numeric derivative w.r.t the next globalJ. Otherwise,
560 // the next "origResidual" will be incorrect.
561 updateCoupledVariables();
562 }
563 }
564};
565
572template<std::size_t id, class TypeTag, class Assembler>
573class SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
574: public SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler,
575 SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, false>, false >
576{
577 using ThisType = SubDomainBoxLocalAssembler<id, TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>;
578 using ParentType = SubDomainBoxLocalAssemblerBase<id, TypeTag, Assembler, ThisType, /*implicit=*/false>;
581 using ElementResidualVector = typename ParentType::LocalResidual::ElementResidualVector;
582 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
584 using FVElementGeometry = typename GridGeometry::LocalView;
585 using GridView = typename GridGeometry::GridView;
586 using Element = typename GridView::template Codim<0>::Entity;
587
589 enum { dim = GridView::dimension };
590
591 static constexpr bool enableGridFluxVarsCache = getPropValue<TypeTag, Properties::EnableGridFluxVariablesCache>();
592 static constexpr bool enableGridVolVarsCache = getPropValue<TypeTag, Properties::EnableGridVolumeVariablesCache>();
593 static constexpr auto domainI = Dune::index_constant<id>();
594
595public:
597
604 template<class JacobianMatrixDiagBlock, class GridVariables>
605 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrixDiagBlock& A, GridVariables& gridVariables)
606 {
607 if (this->assembler().isStationaryProblem())
608 DUNE_THROW(Dune::InvalidStateException, "Using explicit jacobian assembler with stationary local residual");
609
610 // get some aliases for convenience
611 const auto& element = this->element();
612 const auto& fvGeometry = this->fvGeometry();
613 const auto& curSol = this->curSol()[domainI];
614 auto&& curElemVolVars = this->curElemVolVars();
615
616 // get the vector of the actual (undeflected) element residuals
617 const auto origResiduals = this->evalLocalResidual();
618 const auto origStorageResiduals = this->evalLocalStorageResidual();
619
621 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
622 // neighboring elements all derivatives are zero. For the assembled element only the storage //
623 // derivatives are non-zero. //
625
626 // create the element solution
627 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
628
629 // create the vector storing the partial derivatives
630 ElementResidualVector partialDerivs(element.subEntities(dim));
631
632 // calculation of the derivatives
633 for (auto&& scv : scvs(fvGeometry))
634 {
635 // dof index and corresponding actual pri vars
636 const auto dofIdx = scv.dofIndex();
637 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
638 const VolumeVariables origVolVars(curVolVars);
639
640 // calculate derivatives w.r.t to the privars at the dof at hand
641 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
642 {
643 partialDerivs = 0.0;
644
645 auto evalStorage = [&](Scalar priVar)
646 {
647 // auto partialDerivsTmp = partialDerivs;
648 elemSol[scv.localDofIndex()][pvIdx] = priVar;
649 curVolVars.update(elemSol, this->problem(), element, scv);
650 return this->evalLocalStorageResidual();
651 };
652
653 // derive the residuals numerically
654 static const NumericEpsilon<Scalar, numEq> eps_{this->problem().paramGroup()};
655 static const int numDiffMethod = getParamFromGroup<int>(this->problem().paramGroup(), "Assembly.NumericDifferenceMethod");
656 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
657 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
658
659 // update the global stiffness matrix with the current partial derivatives
660 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
661 {
662 // A[i][col][eqIdx][pvIdx] is the rate of change of
663 // the residual of equation 'eqIdx' at dof 'i'
664 // depending on the primary variable 'pvIdx' at dof
665 // 'col'.
666 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
667 }
668
669 // restore the original state of the scv's volume variables
670 curVolVars = origVolVars;
671
672 // restore the original element solution
673 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
674 // TODO additional dof dependencies
675 }
676 }
677
678 // return the undeflected residual
679 return origResiduals;
680 }
681
688 template<std::size_t otherId, class JacobianBlock, class GridVariables>
689 void assembleJacobianCoupling(Dune::index_constant<otherId> domainJ, JacobianBlock& A,
690 const ElementResidualVector& res, GridVariables& gridVariables)
691 {}
692};
693
694} // end namespace Dumux
695
696#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An adapter class for local assemblers using numeric differentiation.
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.
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 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:263
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: fvlocalassemblerbase.hh:68
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:279
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:307
ElementVolumeVariables & prevElemVolVars()
The element volume variables of the provious time step.
Definition: fvlocalassemblerbase.hh:267
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:259
const Assembler & assembler() const
The assembler.
Definition: fvlocalassemblerbase.hh:243
ElementFluxVariablesCache & elemFluxVarsCache()
The element flux variables cache.
Definition: fvlocalassemblerbase.hh:271
LocalResidual & localResidual()
The local residual for the current element.
Definition: fvlocalassemblerbase.hh:275
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:247
const SolutionVector & curSol() const
The current solution.
Definition: fvlocalassemblerbase.hh:255
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:58
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: subdomainboxlocalassembler.hh:231
void bindLocalViews()
Prepares all local views necessary for local assembly.
Definition: subdomainboxlocalassembler.hh:275
static constexpr auto domainId
export the domain id of this sub-domain
Definition: subdomainboxlocalassembler.hh:88
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:205
CouplingManager & couplingManager()
return reference to the coupling manager
Definition: subdomainboxlocalassembler.hh:311
ElementResidualVector evalLocalSourceResidual(const Element &neighbor) const
Evaluates the local source term depending on time discretization scheme.
Definition: subdomainboxlocalassembler.hh:226
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Incorporate Dirichlet boundary conditions.
Definition: subdomainboxlocalassembler.hh:244
void assembleResidual(SubSolutionVector &res)
Assemble the residual vector entries only.
Definition: subdomainboxlocalassembler.hh:182
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:164
const Problem & problem() const
return reference to the underlying problem
Definition: subdomainboxlocalassembler.hh:307
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:114
SubDomainBoxLocalAssemblerBase(const Assembler &assembler, const Element &element, const SolutionVector &curSol, CouplingManager &couplingManager)
Definition: subdomainboxlocalassembler.hh:93
The box scheme multidomain local assembler.
Definition: subdomainboxlocalassembler.hh:330
Box scheme multi domain local assembler using numeric differentiation and implicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:342
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:372
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:460
Box scheme multi domain local assembler using numeric differentiation and explicit time discretizatio...
Definition: subdomainboxlocalassembler.hh:576
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:605
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:689
Declares all properties used in Dumux.