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 "fvbaseproperties.hh"
32
33#include <dune/common/fvector.hh>
34
37
39
40#include <vector>
41
42namespace Opm {
43
50template<class TypeTag>
52{
54
59
60 // the history size of the time discretization in number of steps
61 enum { timeDiscHistorySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
62
63 struct DofStore_ {
64 IntensiveQuantities intensiveQuantities[timeDiscHistorySize];
65 const PrimaryVariables* priVars[timeDiscHistorySize];
66 const IntensiveQuantities *thermodynamicHint[timeDiscHistorySize];
67 };
68 using DofVarsVector = std::vector<DofStore_>;
69 using ExtensiveQuantitiesVector = std::vector<ExtensiveQuantities>;
70
77
79 using Element = typename GridView::template Codim<0>::Entity;
80
81 static const unsigned dimWorld = GridView::dimensionworld;
82 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
83
84 using CoordScalar = typename GridView::ctype;
85 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
86
87 // we don't allow copies of element contexts!
89
90public:
94 explicit FvBaseElementContext(const Simulator& simulator)
96 , stencil_(gridView_, simulator.model().dofMapper() )
97 {
98 // remember the simulator object
100 enableStorageCache_ = Parameters::Get<Parameters::EnableStorageCache>();
101 stashedDofIdx_ = -1;
102 focusDofIdx_ = -1;
103 }
104
105 static void *operator new(size_t size)
106 { return aligned_alloc(alignof(FvBaseElementContext), size); }
107
108 static void operator delete(void *ptr)
109 { aligned_free(ptr); }
110
118 void updateAll(const Element& elem)
119 {
120 asImp_().updateStencil(elem);
121 asImp_().updateAllIntensiveQuantities();
122 asImp_().updateAllExtensiveQuantities();
123 }
124
131 void updateStencil(const Element& elem)
132 {
133 // remember the current element
134 elemPtr_ = &elem;
135
136 // update the stencil. the center gradients are quite expensive to calculate and
137 // most models don't need them, so that we only do this if the model explicitly
138 // enables them
139 stencil_.update(elem);
140
141 // resize the arrays containing the flux and the volume variables
142 dofVars_.resize(stencil_.numDof());
143 extensiveQuantities_.resize(stencil_.numInteriorFaces());
144 }
145
152 void updatePrimaryStencil(const Element& elem)
153 {
154 // remember the current element
155 elemPtr_ = &elem;
156
157 // update the finite element geometry
158 stencil_.updatePrimaryTopology(elem);
159
160 dofVars_.resize(stencil_.numPrimaryDof());
161 }
162
169 void updateStencilTopology(const Element& elem)
170 {
171 // remember the current element
172 elemPtr_ = &elem;
173
174 // update the finite element geometry
175 stencil_.updateTopology(elem);
176 }
177
183 {
184 if (!enableStorageCache_) {
185 // if the storage cache is disabled, we need to calculate the storage term
186 // from scratch, i.e. we need the intensive quantities of all of the history.
187 for (unsigned timeIdx = 0; timeIdx < timeDiscHistorySize; ++ timeIdx)
188 asImp_().updateIntensiveQuantities(timeIdx);
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
203 void updateIntensiveQuantities(unsigned timeIdx)
204 { updateIntensiveQuantities_(timeIdx, numDof(timeIdx)); }
205
212 void updatePrimaryIntensiveQuantities(unsigned timeIdx)
213 { updateIntensiveQuantities_(timeIdx, numPrimaryDof(timeIdx)); }
214
225 void updateIntensiveQuantities(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
226 { asImp_().updateSingleIntQuants_(priVars, dofIdx, timeIdx); }
227
233 { asImp_().updateExtensiveQuantities(/*timeIdx=*/0); }
234
242 void updateExtensiveQuantities(unsigned timeIdx)
243 {
244 gradientCalculator_.prepare(/*context=*/asImp_(), timeIdx);
245
246 for (unsigned fluxIdx = 0; fluxIdx < numInteriorFaces(timeIdx); fluxIdx++) {
247 extensiveQuantities_[fluxIdx].update(/*context=*/asImp_(),
248 /*localIndex=*/fluxIdx,
249 timeIdx);
250 }
251 }
252
260 void setFocusDofIndex(unsigned dofIdx)
261 { focusDofIdx_ = dofIdx; }
262
268 unsigned focusDofIndex() const
269 { return focusDofIdx_; }
270
277 { return this->model().linearizer().getLinearizationType(); }
278
282 const Simulator& simulator() const
283 { return *simulatorPtr_; }
284
288 const Problem& problem() const
289 { return simulatorPtr_->problem(); }
290
294 const Model& model() const
295 { return simulatorPtr_->model(); }
296
300 const GridView& gridView() const
301 { return gridView_; }
302
306 const Element& element() const
307 { return *elemPtr_; }
308
312 size_t numDof(unsigned timeIdx) const
313 { return stencil(timeIdx).numDof(); }
314
318 size_t numPrimaryDof(unsigned timeIdx) const
319 { return stencil(timeIdx).numPrimaryDof(); }
320
325 size_t numInteriorFaces(unsigned timeIdx) const
326 { return stencil(timeIdx).numInteriorFaces(); }
327
332 size_t numBoundaryFaces(unsigned timeIdx) const
333 { return stencil(timeIdx).numBoundaryFaces(); }
334
341 const Stencil& stencil(unsigned) const
342 { return stencil_; }
343
352 decltype(auto) pos(unsigned dofIdx, unsigned) const
353 { return stencil_.subControlVolume(dofIdx).globalPos(); }
354
363 unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
364 { return stencil(timeIdx).globalSpaceIndex(dofIdx); }
365
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#endif
418
419 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
420 }
421
430 const IntensiveQuantities *thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
431 {
432 assert(dofIdx < numDof(timeIdx));
433 return dofVars_[dofIdx].thermodynamicHint[timeIdx];
434 }
438 IntensiveQuantities& intensiveQuantities(unsigned dofIdx, unsigned timeIdx)
439 {
440 assert(dofIdx < numDof(timeIdx));
441 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
442 }
443
452 const PrimaryVariables& primaryVars(unsigned dofIdx, unsigned timeIdx) const
453 {
454 assert(dofIdx < numDof(timeIdx));
455 return *dofVars_[dofIdx].priVars[timeIdx];
456 }
457
465 { return stashedDofIdx_ != -1; }
466
473 int stashedDofIdx() const
474 { return stashedDofIdx_; }
475
481 void stashIntensiveQuantities(unsigned dofIdx)
482 {
483 assert(dofIdx < numDof(/*timeIdx=*/0));
484
485 intensiveQuantitiesStashed_ = dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0];
486 priVarsStashed_ = *dofVars_[dofIdx].priVars[/*timeIdx=*/0];
487 stashedDofIdx_ = static_cast<int>(dofIdx);
488 }
489
495 void restoreIntensiveQuantities(unsigned dofIdx)
496 {
497 dofVars_[dofIdx].priVars[/*timeIdx=*/0] = &priVarsStashed_;
498 dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0] = intensiveQuantitiesStashed_;
499 stashedDofIdx_ = -1;
500 }
501
506 const GradientCalculator& gradientCalculator() const
507 { return gradientCalculator_; }
508
517 const ExtensiveQuantities& extensiveQuantities(unsigned fluxIdx, unsigned) const
518 { return extensiveQuantities_[fluxIdx]; }
519
527 { return enableStorageCache_; }
528
532 void setEnableStorageCache(bool yesno)
533 { enableStorageCache_ = yesno; }
534
535private:
536 Implementation& asImp_()
537 { return *static_cast<Implementation*>(this); }
538
539 const Implementation& asImp_() const
540 { return *static_cast<const Implementation*>(this); }
541
542protected:
548 void updateIntensiveQuantities_(unsigned timeIdx, size_t numDof)
549 {
550 // update the intensive quantities for the whole history
551 const SolutionVector& globalSol = model().solution(timeIdx);
552
553 // update the non-gradient quantities
554 for (unsigned dofIdx = 0; dofIdx < numDof; dofIdx++) {
555 unsigned globalIdx = globalSpaceIndex(dofIdx, timeIdx);
556 const PrimaryVariables& dofSol = globalSol[globalIdx];
557 dofVars_[dofIdx].priVars[timeIdx] = &dofSol;
558
559 dofVars_[dofIdx].thermodynamicHint[timeIdx] =
560 model().thermodynamicHint(globalIdx, timeIdx);
561
562 const auto *cachedIntQuants = model().cachedIntensiveQuantities(globalIdx, timeIdx);
563 if (cachedIntQuants) {
564 dofVars_[dofIdx].intensiveQuantities[timeIdx] = *cachedIntQuants;
565 }
566 else {
567 updateSingleIntQuants_(dofSol, dofIdx, timeIdx);
568 model().updateCachedIntensiveQuantities(dofVars_[dofIdx].intensiveQuantities[timeIdx],
569 globalIdx,
570 timeIdx);
571 }
572 }
573 }
574
575 void updateSingleIntQuants_(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
576 {
577#ifndef NDEBUG
578 if (enableStorageCache_ && timeIdx != 0 && problem().recycleFirstIterationStorage())
579 throw std::logic_error("If caching of the storage term is enabled, only the intensive quantities "
580 "for the most-recent substep (i.e. time index 0) are available!");
581#endif
582
583 dofVars_[dofIdx].priVars[timeIdx] = &priVars;
584 dofVars_[dofIdx].intensiveQuantities[timeIdx].update(/*context=*/asImp_(), dofIdx, timeIdx);
585 }
586
587 IntensiveQuantities intensiveQuantitiesStashed_;
588 PrimaryVariables priVarsStashed_;
589
590 GradientCalculator gradientCalculator_;
591
592 std::vector<DofStore_, aligned_allocator<DofStore_, alignof(DofStore_)> > dofVars_;
593 std::vector<ExtensiveQuantities, aligned_allocator<ExtensiveQuantities, alignof(ExtensiveQuantities)> > extensiveQuantities_;
594
595 const Simulator *simulatorPtr_;
596 const Element *elemPtr_;
597 const GridView gridView_;
598 Stencil stencil_;
599
603};
604
605} // namespace Opm
606
607#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:52
int stashedDofIdx() const
Return the (local) index of the DOF for which the primary variables were stashed.
Definition: fvbaseelementcontext.hh:473
std::vector< ExtensiveQuantities, aligned_allocator< ExtensiveQuantities, alignof(ExtensiveQuantities)> > extensiveQuantities_
Definition: fvbaseelementcontext.hh:593
size_t numPrimaryDof(unsigned timeIdx) const
Return the number of primary degrees of freedom of the current element.
Definition: fvbaseelementcontext.hh:318
bool enableStorageCache() const
Returns true iff the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:526
void updateAllExtensiveQuantities()
Compute the extensive quantities of all sub-control volume faces of the current element for all time ...
Definition: fvbaseelementcontext.hh:232
GradientCalculator gradientCalculator_
Definition: fvbaseelementcontext.hh:590
IntensiveQuantities intensiveQuantitiesStashed_
Definition: fvbaseelementcontext.hh:587
const IntensiveQuantities * thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
Return the thermodynamic hint for a given local index.
Definition: fvbaseelementcontext.hh:430
PrimaryVariables priVarsStashed_
Definition: fvbaseelementcontext.hh:588
void updateAllIntensiveQuantities()
Compute the intensive quantities of all sub-control volumes of the current element for all time indic...
Definition: fvbaseelementcontext.hh:182
const Stencil & stencil(unsigned) const
Return the current finite element geometry.
Definition: fvbaseelementcontext.hh:341
void updateIntensiveQuantities_(unsigned timeIdx, size_t numDof)
Update the first 'n' intensive quantities objects from the primary variables.
Definition: fvbaseelementcontext.hh:548
FvBaseElementContext(const Simulator &simulator)
The constructor.
Definition: fvbaseelementcontext.hh:94
void setEnableStorageCache(bool yesno)
Specifies if the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:532
size_t numDof(unsigned timeIdx) const
Return the number of sub-control volumes of the current element.
Definition: fvbaseelementcontext.hh:312
void updateStencilTopology(const Element &elem)
Update the topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:169
decltype(auto) pos(unsigned dofIdx, unsigned) const
Return the position of a local entities in global coordinates.
Definition: fvbaseelementcontext.hh:352
int focusDofIdx_
Definition: fvbaseelementcontext.hh:601
const Element & element() const
Return the current element.
Definition: fvbaseelementcontext.hh:306
unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
Return the global spatial index for a sub-control volume.
Definition: fvbaseelementcontext.hh:363
LinearizationType linearizationType() const
Returns the linearization type.
Definition: fvbaseelementcontext.hh:276
void updateAll(const Element &elem)
Construct all volume and extensive quantities of an element from scratch.
Definition: fvbaseelementcontext.hh:118
void updatePrimaryIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:212
void updateIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:203
Scalar dofTotalVolume(unsigned dofIdx, unsigned timeIdx) const
Return the total volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:389
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:325
void stashIntensiveQuantities(unsigned dofIdx)
Stash the intensive quantities for a degree of freedom on internal memory.
Definition: fvbaseelementcontext.hh:481
Stencil stencil_
Definition: fvbaseelementcontext.hh:598
bool haveStashedIntensiveQuantities() const
Returns true if no intensive quanties are stashed.
Definition: fvbaseelementcontext.hh:464
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:131
const Simulator * simulatorPtr_
Definition: fvbaseelementcontext.hh:595
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:438
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
const GridView gridView_
Definition: fvbaseelementcontext.hh:597
const PrimaryVariables & primaryVars(unsigned dofIdx, unsigned timeIdx) const
Return the primary variables for a given local index.
Definition: fvbaseelementcontext.hh:452
const GradientCalculator & gradientCalculator() const
Return a reference to the gradient calculation class of the chosen spatial discretization.
Definition: fvbaseelementcontext.hh:506
std::vector< DofStore_, aligned_allocator< DofStore_, alignof(DofStore_)> > dofVars_
Definition: fvbaseelementcontext.hh:592
const Simulator & simulator() const
Return a reference to the simulator.
Definition: fvbaseelementcontext.hh:282
const Element * elemPtr_
Definition: fvbaseelementcontext.hh:596
bool enableStorageCache_
Definition: fvbaseelementcontext.hh:602
const Problem & problem() const
Return a reference to the problem.
Definition: fvbaseelementcontext.hh:288
void updateExtensiveQuantities(unsigned timeIdx)
Compute the extensive quantities of all sub-control volume faces of the current element for a single ...
Definition: fvbaseelementcontext.hh:242
const GridView & gridView() const
Return a reference to the grid view.
Definition: fvbaseelementcontext.hh:300
Scalar dofVolume(unsigned dofIdx, unsigned timeIdx) const
Return the element-local volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:376
size_t numBoundaryFaces(unsigned timeIdx) const
Return the number of boundary faces which need to be considered for the flux apporixmation.
Definition: fvbaseelementcontext.hh:332
const ExtensiveQuantities & extensiveQuantities(unsigned fluxIdx, unsigned) const
Return a reference to the extensive quantities of a sub-control volume face.
Definition: fvbaseelementcontext.hh:517
void updatePrimaryStencil(const Element &elem)
Update the primary topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:152
void updateSingleIntQuants_(const PrimaryVariables &priVars, unsigned dofIdx, unsigned timeIdx)
Definition: fvbaseelementcontext.hh:575
void setFocusDofIndex(unsigned dofIdx)
Sets the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:260
unsigned focusDofIndex() const
Returns the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:268
int stashedDofIdx_
Definition: fvbaseelementcontext.hh:600
void restoreIntensiveQuantities(unsigned dofIdx)
Restores the intensive quantities for a degree of freedom from internal memory.
Definition: fvbaseelementcontext.hh:495
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:225
const Model & model() const
Return a reference to the model.
Definition: fvbaseelementcontext.hh:294
Definition: alignedallocator.hh:114
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:37
void aligned_free(void *ptr) noexcept
Definition: alignedallocator.hh:106
void * aligned_alloc(std::size_t alignment, std::size_t size) noexcept
Definition: alignedallocator.hh:92
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:235
Definition: linearizationtype.hh:35