3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pq1bubblelocalassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_ASSEMBLY_PQ1BUBBLE_LOCAL_ASSEMBLER_HH
26#define DUMUX_ASSEMBLY_PQ1BUBBLE_LOCAL_ASSEMBLER_HH
27
28#include <dune/grid/common/gridenums.hh>
29
40
41#include "volvardeflectionhelper_.hh"
42
43namespace Dumux {
44
46
48{
49 template<class... Args>
50 constexpr void operator()(Args&&...) const {}
51};
52
53template<class X, class Y>
54using Impl = std::conditional_t<!std::is_same_v<X, void>, X, Y>;
55
56} // end namespace Detail
57
67template<class TypeTag, class Assembler, class Implementation, bool implicit>
68class PQ1BubbleLocalAssemblerBase : public FVLocalAssemblerBase<TypeTag, Assembler, Implementation, implicit>
69{
76
77 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
78 static constexpr int dim = GridView::dimension;
79
80public:
81
82 using ParentType::ParentType;
83
85 {
87 this->elemBcTypes().update(this->asImp_().problem(), this->element(), this->fvGeometry());
88 }
89
94 template <class PartialReassembler = DefaultPartialReassembler, class CouplingFunction = Detail::PQ1Bubble::NoOperator>
95 void assembleJacobianAndResidual(JacobianMatrix& jac, SolutionVector& res, GridVariables& gridVariables,
96 const PartialReassembler* partialReassembler,
97 const CouplingFunction& maybeAssembleCouplingBlocks = CouplingFunction())
98 {
99 static_assert(!std::decay_t<decltype(this->asImp_().problem())>::enableInternalDirichletConstraints(), "Internal Dirichlet constraints are currently not implemented for face-centered staggered models!");
100
101 this->asImp_().bindLocalViews();
102 const auto& gridGeometry = this->asImp_().problem().gridGeometry();
103 const auto eIdxGlobal = gridGeometry.elementMapper().index(this->element());
104 if (partialReassembler
105 && partialReassembler->elementColor(eIdxGlobal) == EntityColor::green)
106 {
107 const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
108 for (const auto& scv : scvs(this->fvGeometry()))
109 res[scv.dofIndex()] += residual[scv.localDofIndex()];
110
111 // assemble the coupling blocks for coupled models (does nothing if not coupled)
112 maybeAssembleCouplingBlocks(residual);
113 }
114 else if (!this->elementIsGhost())
115 {
116 const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
117 for (const auto& scv : scvs(this->fvGeometry()))
118 res[scv.dofIndex()] += residual[scv.localDofIndex()];
119
120 // assemble the coupling blocks for coupled models (does nothing if not coupled)
121 maybeAssembleCouplingBlocks(residual);
122 }
123 else
124 {
125 assert(this->elementIsGhost());
126
127 // handle vertex dofs
128 const auto numVerticesLocal = this->element().subEntities(dim);
129 for (int i = 0; i < numVerticesLocal; ++i)
130 {
131 // do not change the non-ghost vertices
132 auto vertex = this->element().template subEntity<dim>(i);
133 if (vertex.partitionType() == Dune::InteriorEntity || vertex.partitionType() == Dune::BorderEntity)
134 continue;
135
136 // set main diagonal entries for the vertex
137 const auto dofIndex = gridGeometry.dofMapper().index(vertex);
138
139 // this might be a vector-valued dof
140 typedef typename JacobianMatrix::block_type BlockType;
141 BlockType &J = jac[dofIndex][dofIndex];
142 for (int j = 0; j < BlockType::rows; ++j)
143 J[j][j] = 1.0;
144
145 // set residual for the ghost vertex dof
146 res[dofIndex] = 0;
147 }
148
149 // handle element dof
150 const auto elemDofIndex = gridGeometry.dofMapper().index(this->element());
151
152 // this might be a vector-valued dof
153 typedef typename JacobianMatrix::block_type BlockType;
154 BlockType &J = jac[elemDofIndex][elemDofIndex];
155 for (int j = 0; j < BlockType::rows; ++j)
156 J[j][j] = 1.0;
157
158 // set residual for the ghost element dof
159 res[elemDofIndex] = 0;
160 }
161
162 auto applyDirichlet = [&] (const auto& scvI,
163 const auto& dirichletValues,
164 const auto eqIdx,
165 const auto pvIdx)
166 {
167 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
168
169 auto& row = jac[scvI.dofIndex()];
170 for (auto col = row.begin(); col != row.end(); ++col)
171 row[col.index()][eqIdx] = 0.0;
172
173 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
174
175
176 // if a periodic dof has Dirichlet values also apply the same Dirichlet values to the other dof TODO periodic
177 // if (this->assembler().gridGeometry().dofOnPeriodicBoundary(scvI.dofIndex()))
178 // {
179 // const auto periodicDof = this->assembler().gridGeometry().periodicallyMappedDof(scvI.dofIndex());
180 // res[periodicDof][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
181 // const auto end = jac[periodicDof].end();
182 // for (auto it = jac[periodicDof].begin(); it != end; ++it)
183 // (*it) = periodicDof != it.index() ? 0.0 : 1.0;
184 // }
185 };
186
187 this->asImp_().enforceDirichletConstraints(applyDirichlet);
188 }
189
194 void assembleJacobian(JacobianMatrix& jac, GridVariables& gridVariables)
195 {
196 this->asImp_().bindLocalViews();
197 this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables); // forward to the internal implementation
198
199 auto applyDirichlet = [&] (const auto& scvI,
200 const auto& dirichletValues,
201 const auto eqIdx,
202 const auto pvIdx)
203 {
204 auto& row = jac[scvI.dofIndex()];
205 for (auto col = row.begin(); col != row.end(); ++col)
206 row[col.index()][eqIdx] = 0.0;
207
208 jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
209 };
210
211 this->asImp_().enforceDirichletConstraints(applyDirichlet);
212 }
213
217 void assembleResidual(SolutionVector& res)
218 {
219 this->asImp_().bindLocalViews();
220 const auto residual = this->evalLocalResidual();
221
222 for (const auto& scv : scvs(this->fvGeometry()))
223 res[scv.dofIndex()] += residual[scv.localDofIndex()];
224
225 auto applyDirichlet = [&] (const auto& scvI,
226 const auto& dirichletValues,
227 const auto eqIdx,
228 const auto pvIdx)
229 {
230 res[scvI.dofIndex()][eqIdx] = this->curElemVolVars()[scvI].priVars()[pvIdx] - dirichletValues[pvIdx];
231 };
232
233 this->asImp_().enforceDirichletConstraints(applyDirichlet);
234 }
235
239 template<typename ApplyFunction>
240 void enforceDirichletConstraints(const ApplyFunction& applyDirichlet)
241 {
242 // enforce Dirichlet boundary conditions
243 this->asImp_().evalDirichletBoundaries(applyDirichlet);
244 // take care of internal Dirichlet constraints (if enabled)
245 this->asImp_().enforceInternalDirichletConstraints(applyDirichlet);
246 }
247
251 template< typename ApplyDirichletFunctionType >
252 void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
253 {
254 // enforce Dirichlet boundaries by overwriting partial derivatives with 1 or 0
255 // and set the residual to (privar - dirichletvalue)
256 if (this->elemBcTypes().hasDirichlet())
257 {
258 for (const auto& scvI : scvs(this->fvGeometry()))
259 {
260 const auto bcTypes = this->elemBcTypes().get(this->fvGeometry(), scvI);
261 if (bcTypes.hasDirichlet())
262 {
263 const auto dirichletValues = this->asImp_().problem().dirichlet(this->element(), scvI);
264
265 // set the Dirichlet conditions in residual and jacobian
266 for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
267 {
268 if (bcTypes.isDirichlet(eqIdx))
269 {
270 const auto pvIdx = bcTypes.eqToDirichletIndex(eqIdx);
271 assert(0 <= pvIdx && pvIdx < numEq);
272 applyDirichlet(scvI, dirichletValues, eqIdx, pvIdx);
273 }
274 }
275 }
276 }
277 }
278 }
279
284 template<class... Args>
285 void maybeUpdateCouplingContext(Args&&...) {}
286
291 template<class... Args>
293};
294
303template<class TypeTag, class Assembler, DiffMethod diffMethod = DiffMethod::numeric, bool implicit = true, class Implementation = void>
305
311template<class TypeTag, class Assembler, class Implementation>
312class PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/true, Implementation>
313: public PQ1BubbleLocalAssemblerBase<TypeTag, Assembler,
314 Detail::PQ1Bubble::Impl<Implementation, PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, true, Implementation>>, true>
315{
319 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
321 using FVElementGeometry = typename GridGeometry::LocalView;
322 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
323 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
327
328 static constexpr auto numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
330
331public:
332
334 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
336
343 template <class PartialReassembler = DefaultPartialReassembler>
344 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
345 const PartialReassembler* partialReassembler = nullptr)
346 {
347 // get some aliases for convenience
348 const auto& problem = this->asImp_().problem();
349 const auto& element = this->element();
350 const auto& fvGeometry = this->fvGeometry();
351 const auto& curSol = this->asImp_().curSol();
352
353 auto&& curElemVolVars = this->curElemVolVars();
354
355 // get the vector of the actual element residuals
356 const auto origResiduals = this->evalLocalResidual();
357
359 // Calculate derivatives of all dofs in the element with respect to the dofs in the stencil. //
361
362 // if all volvars in the stencil have to be updated or if it's enough to only update the
363 // volVars for the scv whose associated dof has been deflected
364 static const bool updateAllVolVars = getParamFromGroup<bool>(
365 this->asImp_().problem().paramGroup(), "Assembly.PQ1BubbleVolVarsDependOnAllElementDofs", false
366 );
367
368 // create the element solution
369 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
370
371 // one residual per element facet
372 const auto numElementResiduals = fvGeometry.numScv();
373
374 // create the vector storing the partial derivatives
375 ElementResidualVector partialDerivs(numElementResiduals);
376
377 Detail::VolVarsDeflectionHelper deflectionHelper(
378 [&] (const auto& scv) -> VolumeVariables& {
379 return this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
380 },
381 fvGeometry,
382 updateAllVolVars
383 );
384
385 // calculation of the derivatives
386 for (const auto& scv : scvs(fvGeometry))
387 {
388 // dof index and corresponding actual pri vars
389 const auto dofIdx = scv.dofIndex();
390 deflectionHelper.setCurrent(scv);
391
392 // calculate derivatives w.r.t to the privars at the dof at hand
393 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
394 {
395 partialDerivs = 0.0;
396
397 auto evalResiduals = [&](Scalar priVar)
398 {
399 // update the volume variables and compute element residual
400 elemSol[scv.localDofIndex()][pvIdx] = priVar;
401 deflectionHelper.deflect(elemSol, scv, problem);
402 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
403 return this->evalLocalResidual();
404 };
405
406 // derive the residuals numerically
407 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
408 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
409 NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
410 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
411
412 // update the global stiffness matrix with the current partial derivatives
413 for (const auto& scvJ : scvs(fvGeometry))
414 {
415 // TODO partial reassembly
416 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
417 {
418 // A[i][col][eqIdx][pvIdx] is the rate of change of
419 // the residual of equation 'eqIdx' at dof 'i'
420 // depending on the primary variable 'pvIdx' at dof
421 // 'col'.
422 A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
423 }
424 }
425
426 // restore the original state of the scv's volume variables
427 deflectionHelper.restore(scv);
428
429 // restore the original element solution
430 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
431 this->asImp_().maybeUpdateCouplingContext(scv, elemSol, pvIdx);
432 }
433 }
434
435 return origResiduals;
436 }
437};
438
439
445template<class TypeTag, class Assembler>
446class PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, /*implicit=*/false>
447: public PQ1BubbleLocalAssemblerBase<TypeTag, Assembler,
448 PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::numeric, false>, false>
449{
453 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
457 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
459
460 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
461
462public:
464
471 template <class PartialReassembler = DefaultPartialReassembler>
472 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
473 const PartialReassembler* partialReassembler = nullptr)
474 {
475 if (partialReassembler)
476 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
477
478 // get some aliases for convenience
479 const auto& problem = this->asImp_().problem();
480 const auto& element = this->element();
481 const auto& fvGeometry = this->fvGeometry();
482 const auto& curSol = this->asImp_().curSol();
483 auto&& curElemVolVars = this->curElemVolVars();
484
485 // get the vecor of the actual element residuals
486 const auto origResiduals = this->evalLocalResidual();
487 const auto origStorageResiduals = this->evalLocalStorageResidual();
488
490 // //
491 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
492 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
493 // actual element. In the actual element we evaluate the derivative of the entire residual. //
494 // //
496
497 // create the element solution
498 auto elemSol = elementSolution(element, curSol, fvGeometry.gridGeometry());
499
500 // create the vector storing the partial derivatives
501 ElementResidualVector partialDerivs(fvGeometry.numScv());
502
503 // calculation of the derivatives
504 for (auto&& scv : scvs(fvGeometry))
505 {
506 // dof index and corresponding actual pri vars
507 const auto dofIdx = scv.dofIndex();
508 auto& curVolVars = this->getVolVarAccess(gridVariables.curGridVolVars(), curElemVolVars, scv);
509 const VolumeVariables origVolVars(curVolVars);
510
511 // calculate derivatives w.r.t to the privars at the dof at hand
512 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
513 {
514 partialDerivs = 0.0;
515
516 auto evalStorage = [&](Scalar priVar)
517 {
518 // auto partialDerivsTmp = partialDerivs;
519 elemSol[scv.localDofIndex()][pvIdx] = priVar;
520 curVolVars.update(elemSol, problem, element, scv);
521 return this->evalLocalStorageResidual();
522 };
523
524 // derive the residuals numerically
525 static const NumericEpsilon<Scalar, numEq> eps_{problem.paramGroup()};
526 static const int numDiffMethod = getParamFromGroup<int>(problem.paramGroup(), "Assembly.NumericDifferenceMethod");
527 NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origStorageResiduals,
528 eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
529
530 // update the global stiffness matrix with the current partial derivatives
531 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
532 {
533 // A[i][col][eqIdx][pvIdx] is the rate of change of
534 // the residual of equation 'eqIdx' at dof 'i'
535 // depending on the primary variable 'pvIdx' at dof
536 // 'col'.
537 A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
538 }
539
540 // restore the original state of the scv's volume variables
541 curVolVars = origVolVars;
542
543 // restore the original element solution
544 elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
545 // TODO additional dof dependencies
546 }
547 }
548 return origResiduals;
549 }
550};
551
557template<class TypeTag, class Assembler>
558class PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/true>
559: public PQ1BubbleLocalAssemblerBase<TypeTag, Assembler,
560 PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, true>, true>
561{
567 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
569
571
572public:
574
581 template <class PartialReassembler = DefaultPartialReassembler>
582 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
583 const PartialReassembler* partialReassembler = nullptr)
584 {
585 if (partialReassembler)
586 DUNE_THROW(Dune::NotImplemented, "partial reassembly for analytic differentiation");
587
588 // get some aliases for convenience
589 const auto& element = this->element();
590 const auto& fvGeometry = this->fvGeometry();
591 const auto& problem = this->asImp_().problem();
592 const auto& curElemVolVars = this->curElemVolVars();
593 const auto& elemFluxVarsCache = this->elemFluxVarsCache();
594
595 // get the vecor of the actual element residuals
596 const auto origResiduals = this->evalLocalResidual();
597
599 // //
600 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
601 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
602 // actual element. In the actual element we evaluate the derivative of the entire residual. //
603 // //
605
606 // calculation of the source and storage derivatives
607 for (const auto& scv : scvs(fvGeometry))
608 {
609 // dof index and corresponding actual pri vars
610 const auto dofIdx = scv.dofIndex();
611 const auto& volVars = curElemVolVars[scv];
612
613 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
614 // only if the problem is instationary we add derivative of storage term
615 // TODO if e.g. porosity depends on all dofs in the element, we would have off-diagonal matrix entries!?
616 if (!this->assembler().isStationaryProblem())
617 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
618 problem,
619 element,
620 fvGeometry,
621 volVars,
622 scv);
623
624 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
625 // add source term derivatives
626 this->localResidual().addSourceDerivatives(A[dofIdx][dofIdx],
627 problem,
628 element,
629 fvGeometry,
630 volVars,
631 scv);
632 }
633
634 // localJacobian[scvIdx][otherScvIdx][eqIdx][priVarIdx] of the fluxes
635 for (const auto& scvf : scvfs(fvGeometry))
636 {
637 if (!scvf.boundary())
638 {
639 // add flux term derivatives
640 this->localResidual().addFluxDerivatives(A,
641 problem,
642 element,
643 fvGeometry,
644 curElemVolVars,
645 elemFluxVarsCache,
646 scvf);
647 }
648
649 // the boundary gets special treatment to simplify
650 // for the user
651 else
652 {
653 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
654 if (this->elemBcTypes().get(fvGeometry, insideScv).hasNeumann())
655 {
656 // add flux term derivatives
657 this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
658 problem,
659 element,
660 fvGeometry,
661 curElemVolVars,
662 elemFluxVarsCache,
663 scvf);
664 }
665 }
666 }
667
668 return origResiduals;
669 }
670};
671
677template<class TypeTag, class Assembler>
678class PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, /*implicit=*/false>
679: public PQ1BubbleLocalAssemblerBase<TypeTag, Assembler,
680 PQ1BubbleLocalAssembler<TypeTag, Assembler, DiffMethod::analytic, false>, false>
681{
687 using ElementResidualVector = typename LocalResidual::ElementResidualVector;
689
690 static constexpr int numEq = GetPropType<TypeTag, Properties::ModelTraits>::numEq();
691
692public:
694
701 template <class PartialReassembler = DefaultPartialReassembler>
702 ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix& A, GridVariables& gridVariables,
703 const PartialReassembler* partialReassembler = nullptr)
704 {
705 if (partialReassembler)
706 DUNE_THROW(Dune::NotImplemented, "partial reassembly for explicit time discretization");
707
708 // get some aliases for convenience
709 const auto& element = this->element();
710 const auto& fvGeometry = this->fvGeometry();
711 const auto& problem = this->asImp_().problem();
712 const auto& curElemVolVars = this->curElemVolVars();
713
714 // get the vector of the actual element residuals
715 const auto origResiduals = this->evalLocalResidual();
716
718 // //
719 // Calculate derivatives of all dofs in stencil with respect to the dofs in the element. In the //
720 // neighboring elements we do so by computing the derivatives of the fluxes which depend on the //
721 // actual element. In the actual element we evaluate the derivative of the entire residual. //
722 // //
724
725 // calculation of the source and storage derivatives
726 for (const auto& scv : scvs(fvGeometry))
727 {
728 // dof index and corresponding actual pri vars
729 const auto dofIdx = scv.dofIndex();
730 const auto& volVars = curElemVolVars[scv];
731
732 // derivative of this scv residual w.r.t the d.o.f. of the same scv (because of mass lumping)
733 // only if the problem is instationary we add derivative of storage term
734 this->localResidual().addStorageDerivatives(A[dofIdx][dofIdx],
735 problem,
736 element,
737 fvGeometry,
738 volVars,
739 scv);
740 }
741
742 return origResiduals;
743 }
744
745}; // explicit PQ1BubbleLocalAssembler with analytic Jacobian
746
747} // end namespace Dumux
748
749#endif
An enum class to define various differentiation methods available in order to compute the derivatives...
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
An adapter class for local assemblers using numeric differentiation.
An assembler for Jacobian and residual contribution per element (box method)
Detects which entries in the Jacobian have to be recomputed.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
DiffMethod
Differentiation methods in order to compute the derivatives of the residual i.e. the entries in the j...
Definition: diffmethod.hh:37
@ green
does not need to be reassembled
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
CVFE< CVFEMethods::PQ1Bubble > PQ1Bubble
Definition: method.hh:97
A base class for all local assemblers.
Definition: fvlocalassemblerbase.hh:48
void bindLocalViews()
Convenience function bind and prepare all relevant variables required for the evaluation of the local...
Definition: fvlocalassemblerbase.hh:185
ElementVolumeVariables & curElemVolVars()
The current element volume variables.
Definition: fvlocalassemblerbase.hh:265
ElementBoundaryTypes & elemBcTypes()
The element's boundary types.
Definition: fvlocalassemblerbase.hh:281
Implementation & asImp_()
Definition: fvlocalassemblerbase.hh:309
ElementResidualVector evalLocalResidual() const
Convenience function to evaluate the complete local residual for the current element....
Definition: fvlocalassemblerbase.hh:120
const Problem & problem() const
The problem.
Definition: fvlocalassemblerbase.hh:241
FVElementGeometry & fvGeometry()
The global finite volume geometry.
Definition: fvlocalassemblerbase.hh:261
bool elementIsGhost() const
Returns if element is a ghost entity.
Definition: fvlocalassemblerbase.hh:253
const Element & element() const
The current element.
Definition: fvlocalassemblerbase.hh:249
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:503
Definition: pq1bubblelocalassembler.hh:48
constexpr void operator()(Args &&...) const
Definition: pq1bubblelocalassembler.hh:50
A base class for all local pq1bubble assemblers.
Definition: pq1bubblelocalassembler.hh:69
void assembleJacobian(JacobianMatrix &jac, GridVariables &gridVariables)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: pq1bubblelocalassembler.hh:194
void assembleJacobianAndResidual(JacobianMatrix &jac, SolutionVector &res, GridVariables &gridVariables, const PartialReassembler *partialReassembler, const CouplingFunction &maybeAssembleCouplingBlocks=CouplingFunction())
Computes the derivatives with respect to the given element and adds them to the global matrix....
Definition: pq1bubblelocalassembler.hh:95
void enforceDirichletConstraints(const ApplyFunction &applyDirichlet)
Enforce Dirichlet constraints.
Definition: pq1bubblelocalassembler.hh:240
void maybeUpdateCouplingContext(Args &&...)
Update the coupling context for coupled models.
Definition: pq1bubblelocalassembler.hh:285
void maybeEvalAdditionalDomainDerivatives(Args &&...)
Update the additional domain derivatives for coupled models.
Definition: pq1bubblelocalassembler.hh:292
void evalDirichletBoundaries(ApplyDirichletFunctionType applyDirichlet)
Evaluates Dirichlet boundaries.
Definition: pq1bubblelocalassembler.hh:252
void bindLocalViews()
Definition: pq1bubblelocalassembler.hh:84
void assembleResidual(SolutionVector &res)
Assemble the residual only.
Definition: pq1bubblelocalassembler.hh:217
An assembler for Jacobian and residual contribution per element (PQ1Bubble methods)
Definition: pq1bubblelocalassembler.hh:304
Face-centered scheme local assembler using numeric differentiation and implicit time discretization.
Definition: pq1bubblelocalassembler.hh:315
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: pq1bubblelocalassembler.hh:344
typename LocalResidual::ElementResidualVector ElementResidualVector
Definition: pq1bubblelocalassembler.hh:334
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: pq1bubblelocalassembler.hh:333
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: pq1bubblelocalassembler.hh:472
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: pq1bubblelocalassembler.hh:582
ElementResidualVector assembleJacobianAndResidualImpl(JacobianMatrix &A, GridVariables &gridVariables, const PartialReassembler *partialReassembler=nullptr)
Computes the derivatives with respect to the given element and adds them to the global matrix.
Definition: pq1bubblelocalassembler.hh:702
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:61
Declares all properties used in Dumux.
The local element solution class for the pq1bubble method.