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
35
36#include <dune/common/fvector.hh>
37
38#include <vector>
39
40namespace Opm {
41
48template<class TypeTag>
50{
52
57
58 // the history size of the time discretization in number of steps
59 enum { timeDiscHistorySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>() };
60
61 struct DofStore_ {
62 IntensiveQuantities intensiveQuantities[timeDiscHistorySize];
63 const PrimaryVariables* priVars[timeDiscHistorySize];
64 const IntensiveQuantities *thermodynamicHint[timeDiscHistorySize];
65 };
66 using DofVarsVector = std::vector<DofStore_>;
67 using ExtensiveQuantitiesVector = std::vector<ExtensiveQuantities>;
68
75
77 using Element = typename GridView::template Codim<0>::Entity;
78
79 static const unsigned dimWorld = GridView::dimensionworld;
80 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
81
82 using CoordScalar = typename GridView::ctype;
83 using GlobalPosition = Dune::FieldVector<CoordScalar, dimWorld>;
84
85 // we don't allow copies of element contexts!
87
88public:
92 explicit FvBaseElementContext(const Simulator& simulator)
94 , stencil_(gridView_, simulator.model().dofMapper() )
95 {
96 // remember the simulator object
98 enableStorageCache_ = Parameters::get<TypeTag, Properties::EnableStorageCache>();
99 stashedDofIdx_ = -1;
100 focusDofIdx_ = -1;
101 }
102
103 static void *operator new(size_t size)
104 { return aligned_alloc(alignof(FvBaseElementContext), size); }
105
106 static void operator delete(void *ptr)
107 { aligned_free(ptr); }
108
116 void updateAll(const Element& elem)
117 {
118 asImp_().updateStencil(elem);
119 asImp_().updateAllIntensiveQuantities();
120 asImp_().updateAllExtensiveQuantities();
121 }
122
129 void updateStencil(const Element& elem)
130 {
131 // remember the current element
132 elemPtr_ = &elem;
133
134 // update the stencil. the center gradients are quite expensive to calculate and
135 // most models don't need them, so that we only do this if the model explicitly
136 // enables them
137 stencil_.update(elem);
138
139 // resize the arrays containing the flux and the volume variables
140 dofVars_.resize(stencil_.numDof());
141 extensiveQuantities_.resize(stencil_.numInteriorFaces());
142 }
143
150 void updatePrimaryStencil(const Element& elem)
151 {
152 // remember the current element
153 elemPtr_ = &elem;
154
155 // update the finite element geometry
156 stencil_.updatePrimaryTopology(elem);
157
158 dofVars_.resize(stencil_.numPrimaryDof());
159 }
160
167 void updateStencilTopology(const Element& elem)
168 {
169 // remember the current element
170 elemPtr_ = &elem;
171
172 // update the finite element geometry
173 stencil_.updateTopology(elem);
174 }
175
181 {
182 if (!enableStorageCache_) {
183 // if the storage cache is disabled, we need to calculate the storage term
184 // from scratch, i.e. we need the intensive quantities of all of the history.
185 for (unsigned timeIdx = 0; timeIdx < timeDiscHistorySize; ++ timeIdx)
186 asImp_().updateIntensiveQuantities(timeIdx);
187 }
188 else
189 // if the storage cache is enabled, we only need to recalculate the storage
190 // term for the most recent point of history (i.e., for the current iterative
191 // solution)
192 asImp_().updateIntensiveQuantities(/*timeIdx=*/0);
193 }
194
201 void updateIntensiveQuantities(unsigned timeIdx)
202 { updateIntensiveQuantities_(timeIdx, numDof(timeIdx)); }
203
210 void updatePrimaryIntensiveQuantities(unsigned timeIdx)
211 { updateIntensiveQuantities_(timeIdx, numPrimaryDof(timeIdx)); }
212
223 void updateIntensiveQuantities(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
224 { asImp_().updateSingleIntQuants_(priVars, dofIdx, timeIdx); }
225
231 { asImp_().updateExtensiveQuantities(/*timeIdx=*/0); }
232
240 void updateExtensiveQuantities(unsigned timeIdx)
241 {
242 gradientCalculator_.prepare(/*context=*/asImp_(), timeIdx);
243
244 for (unsigned fluxIdx = 0; fluxIdx < numInteriorFaces(timeIdx); fluxIdx++) {
245 extensiveQuantities_[fluxIdx].update(/*context=*/asImp_(),
246 /*localIndex=*/fluxIdx,
247 timeIdx);
248 }
249 }
250
258 void setFocusDofIndex(unsigned dofIdx)
259 { focusDofIdx_ = dofIdx; }
260
266 unsigned focusDofIndex() const
267 { return focusDofIdx_; }
268
275 { return this->model().linearizer().getLinearizationType(); }
276
280 const Simulator& simulator() const
281 { return *simulatorPtr_; }
282
286 const Problem& problem() const
287 { return simulatorPtr_->problem(); }
288
292 const Model& model() const
293 { return simulatorPtr_->model(); }
294
298 const GridView& gridView() const
299 { return gridView_; }
300
304 const Element& element() const
305 { return *elemPtr_; }
306
310 size_t numDof(unsigned timeIdx) const
311 { return stencil(timeIdx).numDof(); }
312
316 size_t numPrimaryDof(unsigned timeIdx) const
317 { return stencil(timeIdx).numPrimaryDof(); }
318
323 size_t numInteriorFaces(unsigned timeIdx) const
324 { return stencil(timeIdx).numInteriorFaces(); }
325
330 size_t numBoundaryFaces(unsigned timeIdx) const
331 { return stencil(timeIdx).numBoundaryFaces(); }
332
339 const Stencil& stencil(unsigned) const
340 { return stencil_; }
341
350 decltype(auto) pos(unsigned dofIdx, unsigned) const
351 { return stencil_.subControlVolume(dofIdx).globalPos(); }
352
361 unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
362 { return stencil(timeIdx).globalSpaceIndex(dofIdx); }
363
364
374 Scalar dofVolume(unsigned dofIdx, unsigned timeIdx) const
375 { return stencil(timeIdx).subControlVolume(dofIdx).volume(); }
376
387 Scalar dofTotalVolume(unsigned dofIdx, unsigned timeIdx) const
388 { return model().dofTotalVolume(globalSpaceIndex(dofIdx, timeIdx)); }
389
394 bool onBoundary() const
395 { return element().hasBoundaryIntersections(); }
396
407 const IntensiveQuantities& intensiveQuantities(unsigned dofIdx, unsigned timeIdx) const
408 {
409#ifndef NDEBUG
410 assert(dofIdx < numDof(timeIdx));
411
412 if (enableStorageCache_ && timeIdx != 0 && problem().recycleFirstIterationStorage())
413 throw std::logic_error("If caching of the storage term is enabled, only the intensive quantities "
414 "for the most-recent substep (i.e. time index 0) are available!");
415#endif
416
417 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
418 }
419
428 const IntensiveQuantities *thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
429 {
430 assert(dofIdx < numDof(timeIdx));
431 return dofVars_[dofIdx].thermodynamicHint[timeIdx];
432 }
436 IntensiveQuantities& intensiveQuantities(unsigned dofIdx, unsigned timeIdx)
437 {
438 assert(dofIdx < numDof(timeIdx));
439 return dofVars_[dofIdx].intensiveQuantities[timeIdx];
440 }
441
450 const PrimaryVariables& primaryVars(unsigned dofIdx, unsigned timeIdx) const
451 {
452 assert(dofIdx < numDof(timeIdx));
453 return *dofVars_[dofIdx].priVars[timeIdx];
454 }
455
463 { return stashedDofIdx_ != -1; }
464
471 int stashedDofIdx() const
472 { return stashedDofIdx_; }
473
479 void stashIntensiveQuantities(unsigned dofIdx)
480 {
481 assert(dofIdx < numDof(/*timeIdx=*/0));
482
483 intensiveQuantitiesStashed_ = dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0];
484 priVarsStashed_ = *dofVars_[dofIdx].priVars[/*timeIdx=*/0];
485 stashedDofIdx_ = static_cast<int>(dofIdx);
486 }
487
493 void restoreIntensiveQuantities(unsigned dofIdx)
494 {
495 dofVars_[dofIdx].priVars[/*timeIdx=*/0] = &priVarsStashed_;
496 dofVars_[dofIdx].intensiveQuantities[/*timeIdx=*/0] = intensiveQuantitiesStashed_;
497 stashedDofIdx_ = -1;
498 }
499
504 const GradientCalculator& gradientCalculator() const
505 { return gradientCalculator_; }
506
515 const ExtensiveQuantities& extensiveQuantities(unsigned fluxIdx, unsigned) const
516 { return extensiveQuantities_[fluxIdx]; }
517
525 { return enableStorageCache_; }
526
530 void setEnableStorageCache(bool yesno)
531 { enableStorageCache_ = yesno; }
532
533private:
534 Implementation& asImp_()
535 { return *static_cast<Implementation*>(this); }
536
537 const Implementation& asImp_() const
538 { return *static_cast<const Implementation*>(this); }
539
540protected:
546 void updateIntensiveQuantities_(unsigned timeIdx, size_t numDof)
547 {
548 // update the intensive quantities for the whole history
549 const SolutionVector& globalSol = model().solution(timeIdx);
550
551 // update the non-gradient quantities
552 for (unsigned dofIdx = 0; dofIdx < numDof; dofIdx++) {
553 unsigned globalIdx = globalSpaceIndex(dofIdx, timeIdx);
554 const PrimaryVariables& dofSol = globalSol[globalIdx];
555 dofVars_[dofIdx].priVars[timeIdx] = &dofSol;
556
557 dofVars_[dofIdx].thermodynamicHint[timeIdx] =
558 model().thermodynamicHint(globalIdx, timeIdx);
559
560 const auto *cachedIntQuants = model().cachedIntensiveQuantities(globalIdx, timeIdx);
561 if (cachedIntQuants) {
562 dofVars_[dofIdx].intensiveQuantities[timeIdx] = *cachedIntQuants;
563 }
564 else {
565 updateSingleIntQuants_(dofSol, dofIdx, timeIdx);
566 model().updateCachedIntensiveQuantities(dofVars_[dofIdx].intensiveQuantities[timeIdx],
567 globalIdx,
568 timeIdx);
569 }
570 }
571 }
572
573 void updateSingleIntQuants_(const PrimaryVariables& priVars, unsigned dofIdx, unsigned timeIdx)
574 {
575#ifndef NDEBUG
576 if (enableStorageCache_ && timeIdx != 0 && problem().recycleFirstIterationStorage())
577 throw std::logic_error("If caching of the storage term is enabled, only the intensive quantities "
578 "for the most-recent substep (i.e. time index 0) are available!");
579#endif
580
581 dofVars_[dofIdx].priVars[timeIdx] = &priVars;
582 dofVars_[dofIdx].intensiveQuantities[timeIdx].update(/*context=*/asImp_(), dofIdx, timeIdx);
583 }
584
585 IntensiveQuantities intensiveQuantitiesStashed_;
586 PrimaryVariables priVarsStashed_;
587
588 GradientCalculator gradientCalculator_;
589
590 std::vector<DofStore_, aligned_allocator<DofStore_, alignof(DofStore_)> > dofVars_;
591 std::vector<ExtensiveQuantities, aligned_allocator<ExtensiveQuantities, alignof(ExtensiveQuantities)> > extensiveQuantities_;
592
593 const Simulator *simulatorPtr_;
594 const Element *elemPtr_;
595 const GridView gridView_;
596 Stencil stencil_;
597
601};
602
603} // namespace Opm
604
605#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:50
int stashedDofIdx() const
Return the (local) index of the DOF for which the primary variables were stashed.
Definition: fvbaseelementcontext.hh:471
std::vector< ExtensiveQuantities, aligned_allocator< ExtensiveQuantities, alignof(ExtensiveQuantities)> > extensiveQuantities_
Definition: fvbaseelementcontext.hh:591
size_t numPrimaryDof(unsigned timeIdx) const
Return the number of primary degrees of freedom of the current element.
Definition: fvbaseelementcontext.hh:316
bool enableStorageCache() const
Returns true iff the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:524
void updateAllExtensiveQuantities()
Compute the extensive quantities of all sub-control volume faces of the current element for all time ...
Definition: fvbaseelementcontext.hh:230
GradientCalculator gradientCalculator_
Definition: fvbaseelementcontext.hh:588
IntensiveQuantities intensiveQuantitiesStashed_
Definition: fvbaseelementcontext.hh:585
const IntensiveQuantities * thermodynamicHint(unsigned dofIdx, unsigned timeIdx) const
Return the thermodynamic hint for a given local index.
Definition: fvbaseelementcontext.hh:428
PrimaryVariables priVarsStashed_
Definition: fvbaseelementcontext.hh:586
void updateAllIntensiveQuantities()
Compute the intensive quantities of all sub-control volumes of the current element for all time indic...
Definition: fvbaseelementcontext.hh:180
const Stencil & stencil(unsigned) const
Return the current finite element geometry.
Definition: fvbaseelementcontext.hh:339
void updateIntensiveQuantities_(unsigned timeIdx, size_t numDof)
Update the first 'n' intensive quantities objects from the primary variables.
Definition: fvbaseelementcontext.hh:546
FvBaseElementContext(const Simulator &simulator)
The constructor.
Definition: fvbaseelementcontext.hh:92
void setEnableStorageCache(bool yesno)
Specifies if the cache for the storage term ought to be used for this context.
Definition: fvbaseelementcontext.hh:530
size_t numDof(unsigned timeIdx) const
Return the number of sub-control volumes of the current element.
Definition: fvbaseelementcontext.hh:310
void updateStencilTopology(const Element &elem)
Update the topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:167
decltype(auto) pos(unsigned dofIdx, unsigned) const
Return the position of a local entities in global coordinates.
Definition: fvbaseelementcontext.hh:350
int focusDofIdx_
Definition: fvbaseelementcontext.hh:599
const Element & element() const
Return the current element.
Definition: fvbaseelementcontext.hh:304
unsigned globalSpaceIndex(unsigned dofIdx, unsigned timeIdx) const
Return the global spatial index for a sub-control volume.
Definition: fvbaseelementcontext.hh:361
LinearizationType linearizationType() const
Returns the linearization type.
Definition: fvbaseelementcontext.hh:274
void updateAll(const Element &elem)
Construct all volume and extensive quantities of an element from scratch.
Definition: fvbaseelementcontext.hh:116
void updatePrimaryIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:210
void updateIntensiveQuantities(unsigned timeIdx)
Compute the intensive quantities of all sub-control volumes of the current element for a single time ...
Definition: fvbaseelementcontext.hh:201
Scalar dofTotalVolume(unsigned dofIdx, unsigned timeIdx) const
Return the total volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:387
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:323
void stashIntensiveQuantities(unsigned dofIdx)
Stash the intensive quantities for a degree of freedom on internal memory.
Definition: fvbaseelementcontext.hh:479
Stencil stencil_
Definition: fvbaseelementcontext.hh:596
bool haveStashedIntensiveQuantities() const
Returns true if no intensive quanties are stashed.
Definition: fvbaseelementcontext.hh:462
bool onBoundary() const
Returns whether the current element is on the domain's boundary.
Definition: fvbaseelementcontext.hh:394
void updateStencil(const Element &elem)
Compute the finite volume geometry for an element.
Definition: fvbaseelementcontext.hh:129
const Simulator * simulatorPtr_
Definition: fvbaseelementcontext.hh:593
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:436
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:407
const GridView gridView_
Definition: fvbaseelementcontext.hh:595
const PrimaryVariables & primaryVars(unsigned dofIdx, unsigned timeIdx) const
Return the primary variables for a given local index.
Definition: fvbaseelementcontext.hh:450
const GradientCalculator & gradientCalculator() const
Return a reference to the gradient calculation class of the chosen spatial discretization.
Definition: fvbaseelementcontext.hh:504
std::vector< DofStore_, aligned_allocator< DofStore_, alignof(DofStore_)> > dofVars_
Definition: fvbaseelementcontext.hh:590
const Simulator & simulator() const
Return a reference to the simulator.
Definition: fvbaseelementcontext.hh:280
const Element * elemPtr_
Definition: fvbaseelementcontext.hh:594
bool enableStorageCache_
Definition: fvbaseelementcontext.hh:600
const Problem & problem() const
Return a reference to the problem.
Definition: fvbaseelementcontext.hh:286
void updateExtensiveQuantities(unsigned timeIdx)
Compute the extensive quantities of all sub-control volume faces of the current element for a single ...
Definition: fvbaseelementcontext.hh:240
const GridView & gridView() const
Return a reference to the grid view.
Definition: fvbaseelementcontext.hh:298
Scalar dofVolume(unsigned dofIdx, unsigned timeIdx) const
Return the element-local volume associated with a degree of freedom.
Definition: fvbaseelementcontext.hh:374
size_t numBoundaryFaces(unsigned timeIdx) const
Return the number of boundary faces which need to be considered for the flux apporixmation.
Definition: fvbaseelementcontext.hh:330
const ExtensiveQuantities & extensiveQuantities(unsigned fluxIdx, unsigned) const
Return a reference to the extensive quantities of a sub-control volume face.
Definition: fvbaseelementcontext.hh:515
void updatePrimaryStencil(const Element &elem)
Update the primary topological part of the stencil, but nothing else.
Definition: fvbaseelementcontext.hh:150
void updateSingleIntQuants_(const PrimaryVariables &priVars, unsigned dofIdx, unsigned timeIdx)
Definition: fvbaseelementcontext.hh:573
void setFocusDofIndex(unsigned dofIdx)
Sets the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:258
unsigned focusDofIndex() const
Returns the degree of freedom on which the simulator is currently "focused" on.
Definition: fvbaseelementcontext.hh:266
int stashedDofIdx_
Definition: fvbaseelementcontext.hh:598
void restoreIntensiveQuantities(unsigned dofIdx)
Restores the intensive quantities for a degree of freedom from internal memory.
Definition: fvbaseelementcontext.hh:493
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:223
const Model & model() const
Return a reference to the model.
Definition: fvbaseelementcontext.hh:292
Definition: alignedallocator.hh:114
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:242
Definition: linearizationtype.hh:35