fvbaseelementcontext.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_ELEMENT_CONTEXT_HH
29#define EWOMS_FV_BASE_ELEMENT_CONTEXT_HH
30
31#include <dune/common/fvector.hh>
32
36
38
39#include <array>
40#include <cassert>
41#include <cstddef>
42#include <stdexcept>
43#include <vector>
44
45namespace Opm {
46
53template<class TypeTag>
55{
57
62
63 // the history size of the time discretization in number of steps
64 enum { timeDiscHistorySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
65
66 struct DofStore_ {
67 std::array<IntensiveQuantities, timeDiscHistorySize> intensiveQuantities;
68 std::array<const PrimaryVariables*, timeDiscHistorySize> priVars;
69 std::array<const IntensiveQuantities*, timeDiscHistorySize> thermodynamicHint;
70 };
71 using ExtensiveQuantitiesVector = std::vector<ExtensiveQuantities>;
72
79
81 using Element = typename GridView::template Codim<0>::Entity;
82
83 static const unsigned dimWorld = GridView::dimensionworld;
84
85 using CoordScalar = typename GridView::ctype;
86 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
87
88 // we don't allow copies of element contexts!
90
91public:
95 explicit FvBaseElementContext(const Simulator& simulator)
97 , stencil_(gridView_, simulator.model().dofMapper() )
98 {
99 // remember the simulator object
101 enableStorageCache_ = Parameters::Get<Parameters::EnableStorageCache>();
102 }
103
104 static void *operator new(std::size_t size)
105 { return aligned_alloc(alignof(FvBaseElementContext), size); }
106
107 static void operator delete(void *ptr)
108 { aligned_free(ptr); }
109
117 void updateAll(const Element& elem)
118 {
119 asImp_().updateStencil(elem);
120 asImp_().updateAllIntensiveQuantities();
121 asImp_().updateAllExtensiveQuantities();
122 }
123
130 void updateStencil(const Element& elem)
131 {
132 // remember the current element
133 elemPtr_ = &elem;
134
135 // update the stencil. the center gradients are quite expensive to calculate and
136 // most models don't need them, so that we only do this if the model explicitly
137 // enables them
138 stencil_.update(elem);
139
140 // resize the arrays containing the flux and the volume variables
141 dofVars_.resize(stencil_.numDof());
142 extensiveQuantities_.resize(stencil_.numInteriorFaces());
143 }
144
151 void updatePrimaryStencil(const Element& elem)
152 {
153 // remember the current element
154 elemPtr_ = &elem;
155
156 // update the finite element geometry
157 stencil_.updatePrimaryTopology(elem);
158
159 dofVars_.resize(stencil_.numPrimaryDof());
160 }
161
168 void updateStencilTopology(const Element& elem)
169 {
170 // remember the current element
171 elemPtr_ = &elem;
172
173 // update the finite element geometry
174 stencil_.updateTopology(elem);
175 }
176
182 {
183 if (!enableStorageCache_) {
184 // if the storage cache is disabled, we need to calculate the storage term
185 // from scratch, i.e. we need the intensive quantities of all of the history.
186 for (unsigned timeIdx = 0; timeIdx < timeDiscHistorySize; ++timeIdx) {
187 asImp_().updateIntensiveQuantities(timeIdx);
188 }
189 }
190 else {
191 // if the storage cache is enabled, we only need to recalculate the storage
192 // term for the most recent point of history (i.e., for the current iterative
193 // solution)
194 asImp_().updateIntensiveQuantities(/*timeIdx=*/0);
195 }
196 }
197
204 void updateIntensiveQuantities(unsigned timeIdx)
205 { updateIntensiveQuantities_(timeIdx, numDof(timeIdx)); }
206
213 void updatePrimaryIntensiveQuantities(unsigned timeIdx)
214 { updateIntensiveQuantities_(timeIdx, numPrimaryDof(timeIdx)); }
215
226 void updateIntensiveQuantities(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
227 { asImp_().updateSingleIntQuants_(priVars, dofIdx, timeIdx); }
228
234 { asImp_().updateExtensiveQuantities(/*timeIdx=*/0); }
235
243 void updateExtensiveQuantities(unsigned timeIdx)
244 {
245 gradientCalculator_.prepare(/*context=*/asImp_(), timeIdx);
246
247 for (unsigned fluxIdx = 0; fluxIdx < numInteriorFaces(timeIdx); ++fluxIdx) {
248 extensiveQuantities_[fluxIdx].update(/*context=*/asImp_(),
249 /*localIndex=*/fluxIdx,
250 timeIdx);
251 }
252 }
253
261 void setFocusDofIndex(unsigned dofIdx)
262 { focusDofIdx_ = dofIdx; }
263
269 unsigned focusDofIndex() const
270 { return focusDofIdx_; }
271
278 { return this->model().linearizer().getLinearizationType(); }
279
283 const Simulator& simulator() const
284 { return *simulatorPtr_; }
285
289 const Problem& problem() const
290 { return simulatorPtr_->problem(); }
291
295 const Model& model() const
296 { return simulatorPtr_->model(); }
297
301 const GridView& gridView() const
302 { return gridView_; }
303
307 const Element& element() const
308 { return *elemPtr_; }
309
313 std::size_t numDof(unsigned timeIdx) const
314 { return stencil(timeIdx).numDof(); }
315
319 std::size_t numPrimaryDof(unsigned timeIdx) const
320 { return stencil(timeIdx).numPrimaryDof(); }
321
326 std::size_t numInteriorFaces(unsigned timeIdx) const
327 { return stencil(timeIdx).numInteriorFaces(); }
328
333 std::size_t numBoundaryFaces(unsigned timeIdx) const
334 { return stencil(timeIdx).numBoundaryFaces(); }
335
342 const Stencil& stencil(unsigned) const
343 { return stencil_; }
344
353 decltype(auto) pos(unsigned dofIdx, unsigned) const
354 { return stencil_.subControlVolume(dofIdx).globalPos(); }
355
364 unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
365 { return stencil(timeIdx).globalSpaceIndex(dofIdx); }
366
376 Scalar dofVolume(unsigned dofIdx, unsigned timeIdx) const
377 { return stencil(timeIdx).subControlVolume(dofIdx).volume(); }
378
389 Scalar dofTotalVolume(unsigned dofIdx, unsigned timeIdx) const
390 { return model().dofTotalVolume(globalSpaceIndex(dofIdx, timeIdx)); }
391
396 bool onBoundary() const
397 { return element().hasBoundaryIntersections(); }
398
409 const IntensiveQuantities& intensiveQuantities(unsigned dofIdx, unsigned timeIdx) const
410 {
411#ifndef NDEBUG
412 assert(dofIdx < numDof(timeIdx));
413
414 if (enableStorageCache_ && timeIdx != 0 && problem().recycleFirstIterationStorage()) {
415 throw std::logic_error("If caching of the storage term is enabled, only the intensive quantities "
416 "for the most-recent substep (i.e. time index 0) are available!");
417 }
418#endif
419
420 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
421 }
422
431 const IntensiveQuantities *thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
432 {
433 assert(dofIdx < numDof(timeIdx));
434 return dofVars_[dofIdx].thermodynamicHint[timeIdx];
435 }
436
440 IntensiveQuantities& intensiveQuantities(unsigned dofIdx, unsigned timeIdx)
441 {
442 assert(dofIdx < numDof(timeIdx));
443 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
444 }
445
454 const PrimaryVariables& primaryVars(unsigned dofIdx, unsigned timeIdx) const
455 {
456 assert(dofIdx < numDof(timeIdx));
457 return *dofVars_[dofIdx].priVars[timeIdx];
458 }
459
467 { return stashedDofIdx_ != -1; }
468
475 int stashedDofIdx() const
476 { return stashedDofIdx_; }
477
483 void stashIntensiveQuantities(unsigned dofIdx)
484 {
485 assert(dofIdx < numDof(/*timeIdx=*/0));
486
487 intensiveQuantitiesStashed_ = dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0];
488 priVarsStashed_ = *dofVars_[dofIdx].priVars[/*timeIdx=*/0];
489 stashedDofIdx_ = static_cast<int>(dofIdx);
490 }
491
497 void restoreIntensiveQuantities(unsigned dofIdx)
498 {
499 dofVars_[dofIdx].priVars[/*timeIdx=*/0] = &priVarsStashed_;
500 dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0] = intensiveQuantitiesStashed_;
501 stashedDofIdx_ = -1;
502 }
503
508 const GradientCalculator& gradientCalculator() const
509 { return gradientCalculator_; }
510
519 const ExtensiveQuantities& extensiveQuantities(unsigned fluxIdx, unsigned) const
520 { return extensiveQuantities_[fluxIdx]; }
521
529 { return enableStorageCache_; }
530
534 void setEnableStorageCache(bool yesno)
535 { enableStorageCache_ = yesno; }
536
537private:
538 Implementation& asImp_()
539 { return *static_cast<Implementation*>(this); }
540
541 const Implementation& asImp_() const
542 { return *static_cast<const Implementation*>(this); }
543
544protected:
550 void updateIntensiveQuantities_(unsigned timeIdx, std::size_t numDof)
551 {
552 // update the intensive quantities for the whole history
553 const SolutionVector& globalSol = model().solution(timeIdx);
554
555 // update the non-gradient quantities
556 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
557 unsigned globalIdx = globalSpaceIndex(dofIdx, timeIdx);
558 const PrimaryVariables& dofSol = globalSol[globalIdx];
559 dofVars_[dofIdx].priVars[timeIdx] = &dofSol;
560
561 dofVars_[dofIdx].thermodynamicHint[timeIdx] =
562 model().thermodynamicHint(globalIdx, timeIdx);
563
564 const auto *cachedIntQuants = model().cachedIntensiveQuantities(globalIdx, timeIdx);
565 if (cachedIntQuants) {
566 dofVars_[dofIdx].intensiveQuantities[timeIdx] = *cachedIntQuants;
567 }
568 else {
569 updateSingleIntQuants_(dofSol, dofIdx, timeIdx);
570 model().updateCachedIntensiveQuantities(dofVars_[dofIdx].intensiveQuantities[timeIdx],
571 globalIdx,
572 timeIdx);
573 }
574 }
575 }
576
577 void updateSingleIntQuants_(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
578 {
579#ifndef NDEBUG
580 if (enableStorageCache_ && timeIdx != 0 && problem().recycleFirstIterationStorage()) {
581 throw std::logic_error("If caching of the storage term is enabled, only the intensive quantities "
582 "for the most-recent substep (i.e. time index 0) are available!");
583 }
584#endif
585
586 dofVars_[dofIdx].priVars[timeIdx] = &priVars;
587 dofVars_[dofIdx].intensiveQuantities[timeIdx].update(/*context=*/asImp_(), dofIdx, timeIdx);
588 }
589
590 IntensiveQuantities intensiveQuantitiesStashed_;
591 PrimaryVariables priVarsStashed_;
592
593 GradientCalculator gradientCalculator_;
594
595 template<class T> using AlignedVector = std::vector<T, aligned_allocator<T, alignof(T)>>;
596
599
600 const Simulator* simulatorPtr_{};
601 const Element* elemPtr_{};
602 const GridView gridView_;
603 Stencil stencil_;
604
608};
609
610} // namespace Opm
611
612#endif
This is a stand-alone version of boost::alignment::aligned_allocator from Boost 1....
This class stores an array of IntensiveQuantities objects, one intensive quantities object for each o...
Definition: fvbaseelementcontext.hh:55
int stashedDofIdx() const
Return the (local) index of the DOF for which the primary variables were stashed.
Definition: fvbaseelementcontext.hh:475
std::size_t numBoundaryFaces(unsigned timeIdx) const
Return the number of boundary faces which need to be considered for the flux apporixmation.
Definition: fvbaseelementcontext.hh:333
bool enableStorageCache() const
Returns true iff the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:528
void updateAllExtensiveQuantities()
Compute the extensive quantities of all sub-control volume faces of the current element for all time ...
Definition: fvbaseelementcontext.hh:233
GradientCalculator gradientCalculator_
Definition: fvbaseelementcontext.hh:593
IntensiveQuantities intensiveQuantitiesStashed_
Definition: fvbaseelementcontext.hh:590
std::size_t numDof(unsigned timeIdx) const
Return the number of sub-control volumes of the current element.
Definition: fvbaseelementcontext.hh:313
const IntensiveQuantities * thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
Return the thermodynamic hint for a given local index.
Definition: fvbaseelementcontext.hh:431
PrimaryVariables priVarsStashed_
Definition: fvbaseelementcontext.hh:591
void updateAllIntensiveQuantities()
Compute the intensive quantities of all sub-control volumes of the current element for all time indic...
Definition: fvbaseelementcontext.hh:181
const Stencil & stencil(unsigned) const
Return the current finite element geometry.
Definition: fvbaseelementcontext.hh:342
std::vector< T, aligned_allocator< T, alignof(T)> > AlignedVector
Definition: fvbaseelementcontext.hh:595
FvBaseElementContext(const Simulator &simulator)
The constructor.
Definition: fvbaseelementcontext.hh:95
void setEnableStorageCache(bool yesno)
Specifies if the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:534
void updateStencilTopology(const Element &elem)
Update the topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:168
decltype(auto) pos(unsigned dofIdx, unsigned) const
Return the position of a local entities in global coordinates.
Definition: fvbaseelementcontext.hh:353
int focusDofIdx_
Definition: fvbaseelementcontext.hh:606
const Element & element() const
Return the current element.
Definition: fvbaseelementcontext.hh:307
unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
Return the global spatial index for a sub-control volume.
Definition: fvbaseelementcontext.hh:364
LinearizationType linearizationType() const
Returns the linearization type.
Definition: fvbaseelementcontext.hh:277
void updateAll(const Element &elem)
Construct all volume and extensive quantities of an element from scratch.
Definition: fvbaseelementcontext.hh:117
void updatePrimaryIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:213
void updateIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:204
Scalar dofTotalVolume(unsigned dofIdx, unsigned timeIdx) const
Return the total volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:389
void stashIntensiveQuantities(unsigned dofIdx)
Stash the intensive quantities for a degree of freedom on internal memory.
Definition: fvbaseelementcontext.hh:483
Stencil stencil_
Definition: fvbaseelementcontext.hh:603
bool haveStashedIntensiveQuantities() const
Returns true if no intensive quanties are stashed.
Definition: fvbaseelementcontext.hh:466
bool onBoundary() const
Returns whether the current element is on the domain's boundary.
Definition: fvbaseelementcontext.hh:396
void updateStencil(const Element &elem)
Compute the finite volume geometry for an element.
Definition: fvbaseelementcontext.hh:130
const Simulator * simulatorPtr_
Definition: fvbaseelementcontext.hh:600
IntensiveQuantities & intensiveQuantities(unsigned dofIdx, unsigned timeIdx)
Return a reference to the intensive quantities of a sub-control volume at a given time.
Definition: fvbaseelementcontext.hh:440
const IntensiveQuantities & intensiveQuantities(unsigned dofIdx, unsigned timeIdx) const
Return a reference to the intensive quantities of a sub-control volume at a given time.
Definition: fvbaseelementcontext.hh:409
void updateIntensiveQuantities_(unsigned timeIdx, std::size_t numDof)
Update the first 'n' intensive quantities objects from the primary variables.
Definition: fvbaseelementcontext.hh:550
const GridView gridView_
Definition: fvbaseelementcontext.hh:602
const PrimaryVariables & primaryVars(unsigned dofIdx, unsigned timeIdx) const
Return the primary variables for a given local index.
Definition: fvbaseelementcontext.hh:454
const GradientCalculator & gradientCalculator() const
Return a reference to the gradient calculation class of the chosen spatial discretization.
Definition: fvbaseelementcontext.hh:508
AlignedVector< ExtensiveQuantities > extensiveQuantities_
Definition: fvbaseelementcontext.hh:598
const Simulator & simulator() const
Return a reference to the simulator.
Definition: fvbaseelementcontext.hh:283
std::size_t numInteriorFaces(unsigned timeIdx) const
Return the number of non-boundary faces which need to be considered for the flux apporixmation.
Definition: fvbaseelementcontext.hh:326
AlignedVector< DofStore_ > dofVars_
Definition: fvbaseelementcontext.hh:597
const Element * elemPtr_
Definition: fvbaseelementcontext.hh:601
bool enableStorageCache_
Definition: fvbaseelementcontext.hh:607
std::size_t numPrimaryDof(unsigned timeIdx) const
Return the number of primary degrees of freedom of the current element.
Definition: fvbaseelementcontext.hh:319
const Problem & problem() const
Return a reference to the problem.
Definition: fvbaseelementcontext.hh:289
void updateExtensiveQuantities(unsigned timeIdx)
Compute the extensive quantities of all sub-control volume faces of the current element for a single ...
Definition: fvbaseelementcontext.hh:243
const GridView & gridView() const
Return a reference to the grid view.
Definition: fvbaseelementcontext.hh:301
Scalar dofVolume(unsigned dofIdx, unsigned timeIdx) const
Return the element-local volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:376
const ExtensiveQuantities & extensiveQuantities(unsigned fluxIdx, unsigned) const
Return a reference to the extensive quantities of a sub-control volume face.
Definition: fvbaseelementcontext.hh:519
void updatePrimaryStencil(const Element &elem)
Update the primary topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:151
void updateSingleIntQuants_(const PrimaryVariables &priVars, unsigned dofIdx, unsigned timeIdx)
Definition: fvbaseelementcontext.hh:577
void setFocusDofIndex(unsigned dofIdx)
Sets the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:261
unsigned focusDofIndex() const
Returns the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:269
int stashedDofIdx_
Definition: fvbaseelementcontext.hh:605
void restoreIntensiveQuantities(unsigned dofIdx)
Restores the intensive quantities for a degree of freedom from internal memory.
Definition: fvbaseelementcontext.hh:497
void updateIntensiveQuantities(const PrimaryVariables &priVars, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities of a single sub-control volume of the current element for a single t...
Definition: fvbaseelementcontext.hh:226
const Model & model() const
Return a reference to the model.
Definition: fvbaseelementcontext.hh:295
Definition: alignedallocator.hh:97
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilboundaryratevector.hh:39
void aligned_free(void *ptr) noexcept
Definition: alignedallocator.hh:89
void * aligned_alloc(std::size_t alignment, std::size_t size) noexcept
Definition: alignedallocator.hh:75
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
Definition: linearizationtype.hh:34