fvbaselocalresidual.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_FV_BASE_LOCAL_RESIDUAL_HH
29#define EWOMS_FV_BASE_LOCAL_RESIDUAL_HH
30
31#include <dune/common/classname.hh>
32#include <dune/common/fvector.hh>
33
34#include <dune/grid/common/geometry.hh>
35
36#include <dune/istl/bvector.hh>
37
38#include <opm/material/common/MathToolbox.hpp>
39#include <opm/material/common/Valgrind.hpp>
40
44
45#include <cassert>
46#include <cmath>
47#include <cstddef>
48#include <stdexcept>
49#include <type_traits>
50
51namespace Opm {
52
61template<class TypeTag>
63{
64private:
66
68 using Element = typename GridView::template Codim<0>::Entity;
69
79
80 static constexpr bool useVolumetricResidual = getPropValue<TypeTag, Properties::UseVolumetricResidual>();
81
82 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
83 enum { extensiveStorageTerm = getPropValue<TypeTag, Properties::ExtensiveStorageTerm>() };
84
85 using Toolbox = MathToolbox<Evaluation>;
86 using EvalVector = Dune::FieldVector<Evaluation, numEq>;
87
88public:
89 using LocalEvalBlockVector = Dune::BlockVector<EvalVector, aligned_allocator<EvalVector, alignof(EvalVector)>>;
90
92
93 // copying the local residual class is not a good idea
95
99 static void registerParameters()
100 {}
101
107 { return internalResidual_; }
108
115 const EvalVector& residual(unsigned dofIdx) const
116 { return internalResidual_[dofIdx]; }
117
128 void eval(const Problem& problem, const Element& element)
129 {
130 ElementContext elemCtx(problem);
131 elemCtx.updateAll(element);
132 eval(elemCtx);
133 }
134
144 void eval(ElementContext& elemCtx)
145 {
146 std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
147 internalResidual_.resize(numDof);
148 asImp_().eval(internalResidual_, elemCtx);
149 }
150
159 ElementContext& elemCtx) const
160 {
161 assert(residual.size() == elemCtx.numDof(/*timeIdx=*/0));
162
163 residual = 0.0;
164
165 // evaluate the flux terms
166 asImp_().evalFluxes(residual, elemCtx, /*timeIdx=*/0);
167
168 // evaluate the storage and the source terms
169 asImp_().evalVolumeTerms_(residual, elemCtx);
170
171 // evaluate the boundary conditions
172 asImp_().evalBoundary_(residual, elemCtx, /*timeIdx=*/0);
173
174 if (useVolumetricResidual) {
175 // make the residual volume specific (i.e., make it incorrect mass per cubic
176 // meter instead of total mass)
177 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
178 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
179 if (elemCtx.dofTotalVolume(dofIdx, /*timeIdx=*/0) > 0.0) {
180 // interior DOF
181 const Scalar dofVolume = elemCtx.dofTotalVolume(dofIdx, /*timeIdx=*/0);
182
183 assert(std::isfinite(dofVolume));
184 Valgrind::CheckDefined(dofVolume);
185
186 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
187 residual[dofIdx][eqIdx] /= dofVolume;
188 }
189 }
190 }
191 }
192 }
193
206 const ElementContext& elemCtx,
207 unsigned timeIdx) const
208 {
209 // the derivative of the storage term depends on the current primary variables;
210 // for time indices != 0, the storage term is constant (because these solutions
211 // are not changed by the Newton method!)
212 if (timeIdx == 0) {
213 // calculate the amount of conservation each quantity inside
214 // all primary sub control volumes
215 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
216 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
217 storage[dofIdx] = 0.0;
218
219 // the volume of the associated DOF
220 const Scalar alpha =
221 elemCtx.stencil(timeIdx).subControlVolume(dofIdx).volume() *
222 elemCtx.intensiveQuantities(dofIdx, timeIdx).extrusionFactor();
223
224 // If the degree of freedom which we currently look at is the one at the
225 // center of attention, we need to consider the derivatives for the
226 // storage term, else the storage term is constant w.r.t. the primary
227 // variables of the focused DOF.
228 if (dofIdx == elemCtx.focusDofIndex()) {
229 asImp_().computeStorage(storage[dofIdx],
230 elemCtx,
231 dofIdx,
232 timeIdx);
233
234 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
235 storage[dofIdx][eqIdx] *= alpha;
236 }
237 }
238 else {
239 Dune::FieldVector<Scalar, numEq> tmp;
240 asImp_().computeStorage(tmp,
241 elemCtx,
242 dofIdx,
243 timeIdx);
244
245 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
246 storage[dofIdx][eqIdx] = tmp[eqIdx] * alpha;
247 }
248 }
249 }
250 }
251 else {
252 // for all previous solutions, the storage term does _not_ depend on the
253 // current primary variables, so we use scalars to store it.
254 if (elemCtx.enableStorageCache()) {
255 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
256 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
257 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
258 const auto& cachedStorage = elemCtx.model().cachedStorage(globalDofIdx, timeIdx);
259 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
260 storage[dofIdx][eqIdx] = cachedStorage[eqIdx];
261 }
262 }
263 }
264 else {
265 // calculate the amount of conservation each quantity inside
266 // all primary sub control volumes
267 Dune::FieldVector<Scalar, numEq> tmp;
268 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
269 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
270 tmp = 0.0;
271 asImp_().computeStorage(tmp,
272 elemCtx,
273 dofIdx,
274 timeIdx);
275 tmp *= elemCtx.stencil(timeIdx).subControlVolume(dofIdx).volume() *
276 elemCtx.intensiveQuantities(dofIdx, timeIdx).extrusionFactor();
277
278 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx)
279 storage[dofIdx][eqIdx] = tmp[eqIdx];
280 }
281 }
282 }
283
284#ifndef NDEBUG
285 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
286 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
287 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
288 Valgrind::CheckDefined(storage[dofIdx][eqIdx]);
289 assert(isfinite(storage[dofIdx][eqIdx]));
290 }
291 }
292#endif
293 }
294
303 const ElementContext& elemCtx,
304 unsigned timeIdx) const
305 {
306 RateVector flux;
307
308 const auto& stencil = elemCtx.stencil(timeIdx);
309 // calculate the mass flux over the sub-control volume faces
310 const std::size_t numInteriorFaces = elemCtx.numInteriorFaces(timeIdx);
311 for (unsigned scvfIdx = 0; scvfIdx < numInteriorFaces; ++scvfIdx) {
312 const auto& face = stencil.interiorFace(scvfIdx);
313 const unsigned i = face.interiorIndex();
314 const unsigned j = face.exteriorIndex();
315
316 Valgrind::SetUndefined(flux);
317 asImp_().computeFlux(flux, /*context=*/elemCtx, scvfIdx, timeIdx);
318 Valgrind::CheckDefined(flux);
319#ifndef NDEBUG
320 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
321 assert(isfinite(flux[eqIdx]));
322 }
323#endif
324
325 const Scalar alpha = elemCtx.extensiveQuantities(scvfIdx, timeIdx).extrusionFactor() * face.area();
326 Valgrind::CheckDefined(alpha);
327 assert(alpha > 0.0);
328 assert(isfinite(alpha));
329
330 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
331 flux[eqIdx] *= alpha;
332 }
333
334 // The balance equation for a finite volume is given by
335 //
336 // dStorage/dt + Flux = Source
337 //
338 // where the 'Flux' and the 'Source' terms represent the
339 // mass per second which leaves the finite
340 // volume. Re-arranging this, we get
341 //
342 // dStorage/dt + Flux - Source = 0
343 //
344 // Since the mass flux as calculated by computeFlux() goes out of sub-control
345 // volume i and into sub-control volume j, we need to add the flux to finite
346 // volume i and subtract it from finite volume j
347 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
348 assert(isfinite(flux[eqIdx]));
349 residual[i][eqIdx] += flux[eqIdx];
350 residual[j][eqIdx] -= flux[eqIdx];
351 }
352 }
353
354#if !defined NDEBUG
355 // in debug mode, ensure that the residual is well-defined
356 const std::size_t numDof = elemCtx.numDof(timeIdx);
357 for (unsigned i = 0; i < numDof; ++i) {
358 for (unsigned j = 0; j < numEq; ++j) {
359 assert(isfinite(residual[i][j]));
360 Valgrind::CheckDefined(residual[i][j]);
361 }
362 }
363#endif
364 }
365
367 // The following methods _must_ be overloaded by the actual flow
368 // models!
370
378 void computeStorage(EqVector&,
379 const ElementContext&,
380 unsigned,
381 unsigned) const
382 {
383 throw std::logic_error("Not implemented: The local residual " + Dune::className<Implementation>() +
384 " does not implement the required method 'computeStorage()'");
385 }
386
394 void computeFlux(RateVector&,
395 const ElementContext&,
396 unsigned,
397 unsigned) const
398 {
399 throw std::logic_error("Not implemented: The local residual " + Dune::className<Implementation>() +
400 " does not implement the required method 'computeFlux()'");
401 }
402
409 void computeSource(RateVector&,
410 const ElementContext&,
411 unsigned,
412 unsigned) const
413 {
414 throw std::logic_error("Not implemented: The local residual " + Dune::className<Implementation>() +
415 " does not implement the required method 'computeSource()'");
416 }
417
418protected:
423 const ElementContext& elemCtx,
424 unsigned timeIdx) const
425 {
426 if (!elemCtx.onBoundary()) {
427 return;
428 }
429
430 BoundaryContext boundaryCtx(elemCtx);
431 // move the iterator to the first boundary
432 if (boundaryCtx.intersection(0).neighbor()) {
433 boundaryCtx.increment();
434 }
435
436 // evaluate the boundary for all boundary faces of the current context
437 const std::size_t numBoundaryFaces = boundaryCtx.numBoundaryFaces(/*timeIdx=*/0);
438 for (unsigned faceIdx = 0; faceIdx < numBoundaryFaces; ++faceIdx, boundaryCtx.increment()) {
439 // add the residual of all vertices of the boundary
440 // segment
442 boundaryCtx,
443 faceIdx,
444 timeIdx);
445 }
446
447#if !defined NDEBUG
448 // in debug mode, ensure that the residual and the storage terms are well-defined
449 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
450 for (unsigned i = 0; i < numDof; ++i) {
451 for (unsigned j = 0; j < numEq; ++j) {
452 assert(isfinite(residual[i][j]));
453 Valgrind::CheckDefined(residual[i][j]);
454 }
455 }
456#endif
457 }
458
464 const BoundaryContext& boundaryCtx,
465 unsigned boundaryFaceIdx,
466 unsigned timeIdx) const
467 {
468 BoundaryRateVector values;
469
470 Valgrind::SetUndefined(values);
471 boundaryCtx.problem().boundary(values, boundaryCtx, boundaryFaceIdx, timeIdx);
472 Valgrind::CheckDefined(values);
473
474 const auto& stencil = boundaryCtx.stencil(timeIdx);
475 const unsigned dofIdx = stencil.boundaryFace(boundaryFaceIdx).interiorIndex();
476 const auto& insideIntQuants = boundaryCtx.elementContext().intensiveQuantities(dofIdx, timeIdx);
477 for (unsigned eqIdx = 0; eqIdx < values.size(); ++eqIdx) {
478 values[eqIdx] *= stencil.boundaryFace(boundaryFaceIdx).area() *
479 insideIntQuants.extrusionFactor();
480
481 Valgrind::CheckDefined(values[eqIdx]);
482 assert(isfinite(values[eqIdx]));
483 }
484
485 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
486 residual[dofIdx][eqIdx] += values[eqIdx];
487 }
488 }
489
496 ElementContext& elemCtx) const
497 {
498 EvalVector tmp;
499 EqVector tmp2;
500 RateVector sourceRate;
501
502 tmp = 0.0;
503 tmp2 = 0.0;
504
505 // evaluate the volumetric terms (storage + source terms)
506 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
507 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
508 const Scalar extrusionFactor =
509 elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0).extrusionFactor();
510 Valgrind::CheckDefined(extrusionFactor);
511 assert(isfinite(extrusionFactor));
512 assert(extrusionFactor > 0.0);
513 const Scalar scvVolume =
514 elemCtx.stencil(/*timeIdx=*/0).subControlVolume(dofIdx).volume() * extrusionFactor;
515 Valgrind::CheckDefined(scvVolume);
516 assert(isfinite(scvVolume));
517 assert(scvVolume > 0.0);
518
519 // if the model uses extensive quantities in its storage term, and we use
520 // automatic differention and current DOF is also not the one we currently
521 // focus on, the storage term does not need any derivatives!
522 if (!extensiveStorageTerm &&
523 !std::is_same_v<Scalar, Evaluation> &&
524 dofIdx != elemCtx.focusDofIndex())
525 {
526 asImp_().computeStorage(tmp2, elemCtx, dofIdx, /*timeIdx=*/0);
527 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
528 tmp[eqIdx] = tmp2[eqIdx];
529 }
530 }
531 else {
532 asImp_().computeStorage(tmp, elemCtx, dofIdx, /*timeIdx=*/0);
533 }
534
535#ifndef NDEBUG
536 Valgrind::CheckDefined(tmp);
537 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
538 assert(isfinite(tmp[eqIdx]));
539 }
540#endif
541
542 if (elemCtx.enableStorageCache()) {
543 const auto& model = elemCtx.model();
544 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
545 if (model.newtonMethod().numIterations() == 0 &&
546 !elemCtx.haveStashedIntensiveQuantities())
547 {
548 if (!elemCtx.problem().recycleFirstIterationStorage()) {
549 // we re-calculate the storage term for the solution of the
550 // previous time step from scratch instead of using the one of
551 // the first iteration of the current time step.
552 tmp2 = 0.0;
553 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/1);
554 asImp_().computeStorage(tmp2, elemCtx, dofIdx, /*timeIdx=*/1);
555 }
556 else {
557 // if the storage term is cached and we're in the first iteration
558 // of the time step, use the storage term of the first iteration
559 // as the one as the solution of the last time step (this assumes
560 // that the initial guess for the solution at the end of the time
561 // step is the same as the solution at the beginning of the time
562 // step. This is usually true, but some fancy preprocessing
563 // scheme might invalidate that assumption.)
564 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
565 tmp2[eqIdx] = Toolbox::value(tmp[eqIdx]);
566 }
567 }
568
569 Valgrind::CheckDefined(tmp2);
570
571 model.updateCachedStorage(globalDofIdx, /*timeIdx=*/1, tmp2);
572 }
573 else {
574 // if the mass storage at the beginning of the time step is not cached,
575 // if the storage term is cached and we're not looking at the first
576 // iteration of the time step, we take the cached data.
577 tmp2 = model.cachedStorage(globalDofIdx, /*timeIdx=*/1);
578 Valgrind::CheckDefined(tmp2);
579 }
580 }
581 else {
582 // if the mass storage at the beginning of the time step is not cached,
583 // we re-calculate it from scratch.
584 tmp2 = 0.0;
585 asImp_().computeStorage(tmp2, elemCtx, dofIdx, /*timeIdx=*/1);
586 Valgrind::CheckDefined(tmp2);
587 }
588
589 // Use the implicit Euler time discretization
590 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
591 const double dt = elemCtx.simulator().timeStepSize();
592 assert(dt > 0);
593 tmp[eqIdx] -= tmp2[eqIdx];
594 tmp[eqIdx] *= scvVolume / dt;
595
596 residual[dofIdx][eqIdx] += tmp[eqIdx];
597 }
598
599 Valgrind::CheckDefined(residual[dofIdx]);
600
601 // deal with the source term
602 asImp_().computeSource(sourceRate, elemCtx, dofIdx, /*timeIdx=*/0);
603
604 // if the model uses extensive quantities in its storage term, and we use
605 // automatic differention and current DOF is also not the one we currently
606 // focus on, the storage term does not need any derivatives!
607 if (!extensiveStorageTerm &&
608 !std::is_same_v<Scalar, Evaluation> &&
609 dofIdx != elemCtx.focusDofIndex())
610 {
611 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
612 residual[dofIdx][eqIdx] -= scalarValue(sourceRate[eqIdx])*scvVolume;
613 }
614 }
615 else {
616 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
617 sourceRate[eqIdx] *= scvVolume;
618 residual[dofIdx][eqIdx] -= sourceRate[eqIdx];
619 }
620 }
621
622 Valgrind::CheckDefined(residual[dofIdx]);
623 }
624
625#if !defined NDEBUG
626 // in debug mode, ensure that the residual is well-defined
627 std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
628 for (unsigned i = 0; i < numDof; ++i) {
629 for (unsigned j = 0; j < numEq; ++j) {
630 assert(isfinite(residual[i][j]));
631 Valgrind::CheckDefined(residual[i][j]);
632 }
633 }
634#endif
635 }
636
637private:
638 Implementation& asImp_()
639 { return *static_cast<Implementation*>(this); }
640
641 const Implementation& asImp_() const
642 { return *static_cast<const Implementation*>(this); }
643
644 LocalEvalBlockVector internalResidual_;
645};
646
647} // namespace Opm
648
649#endif
This is a stand-alone version of boost::alignment::aligned_allocator from Boost 1....
Element-wise caculation of the residual matrix for models based on a finite volume spatial discretiza...
Definition: fvbaselocalresidual.hh:63
Dune::BlockVector< EvalVector, aligned_allocator< EvalVector, alignof(EvalVector)> > LocalEvalBlockVector
Definition: fvbaselocalresidual.hh:89
void evalVolumeTerms_(LocalEvalBlockVector &residual, ElementContext &elemCtx) const
Add the change in the storage terms and the source term to the local residual of all sub-control volu...
Definition: fvbaselocalresidual.hh:495
static void registerParameters()
Register all run-time parameters for the local residual.
Definition: fvbaselocalresidual.hh:99
void computeSource(RateVector &, const ElementContext &, unsigned, unsigned) const
Calculate the source term of the equation.
Definition: fvbaselocalresidual.hh:409
FvBaseLocalResidual(const FvBaseLocalResidual &)=delete
void evalFluxes(LocalEvalBlockVector &residual, const ElementContext &elemCtx, unsigned timeIdx) const
Add the flux term to a local residual.
Definition: fvbaselocalresidual.hh:302
void computeStorage(EqVector &, const ElementContext &, unsigned, unsigned) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: fvbaselocalresidual.hh:378
void evalStorage(LocalEvalBlockVector &storage, const ElementContext &elemCtx, unsigned timeIdx) const
Calculate the amount of all conservation quantities stored in all element's sub-control volumes for a...
Definition: fvbaselocalresidual.hh:205
void eval(const Problem &problem, const Element &element)
Compute the local residual, i.e. the deviation of the conservation equations from zero and store the ...
Definition: fvbaselocalresidual.hh:128
void evalBoundarySegment_(LocalEvalBlockVector &residual, const BoundaryContext &boundaryCtx, unsigned boundaryFaceIdx, unsigned timeIdx) const
Evaluate all boundary conditions for a single sub-control volume face to the local residual.
Definition: fvbaselocalresidual.hh:463
void eval(LocalEvalBlockVector &residual, ElementContext &elemCtx) const
Compute the local residual, i.e. the deviation of the conservation equations from zero.
Definition: fvbaselocalresidual.hh:158
void computeFlux(RateVector &, const ElementContext &, unsigned, unsigned) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: fvbaselocalresidual.hh:394
void evalBoundary_(LocalEvalBlockVector &residual, const ElementContext &elemCtx, unsigned timeIdx) const
Evaluate the boundary conditions of an element.
Definition: fvbaselocalresidual.hh:422
const LocalEvalBlockVector & residual() const
Return the result of the eval() call using internal storage.
Definition: fvbaselocalresidual.hh:106
void eval(ElementContext &elemCtx)
Compute the local residual, i.e. the deviation of the conservation equations from zero and store the ...
Definition: fvbaselocalresidual.hh:144
const EvalVector & residual(unsigned dofIdx) const
Return the result of the eval() call using internal storage.
Definition: fvbaselocalresidual.hh:115
Definition: alignedallocator.hh:97
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:39
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:233
This file provides the infrastructure to retrieve run-time parameters.