3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2p2c/sequential/celldata.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_ELEMENTDATA2P2C_HH
26#define DUMUX_ELEMENTDATA2P2C_HH
27
28#include "fluxdata.hh"
29
30namespace Dumux {
42template<class TypeTag>
44{
45private:
49
51
52 enum
53 {
54 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
55 wCompIdx = Indices::wPhaseIdx, nCompIdx = Indices::nPhaseIdx,
56 contiWEqIdx = Indices::contiWEqIdx
57 };
58 enum
59 {
60 numPhases = getPropValue<TypeTag, Properties::NumPhases>(),
61 numComponents = getPropValue<TypeTag, Properties::NumComponents>()
62 };
63protected:
64 // primary variable (phase pressure has to be stored in fluidstate)
65 Scalar massConcentration_[numComponents];
66
67 Scalar mobility_[numPhases];
68
69 //numerical quantities
72 Scalar dv_dp_;
73 Scalar dv_[numComponents];
76
78 Scalar perimeter_;
79
80 std::shared_ptr<FluidState> fluidState_;
82public:
83
86 {
87 for (int i = 0; i < numPhases;i++)
88 {
89 mobility_[i] = 0.0;
90 dv_[i] = 0.0;
91 }
92 volumeError_ = 0.;
94 dv_dp_ = 0.;
96 wasRefined_ = false;
97 globalIdx_ = 0;
98 perimeter_ = 0.;
99 }
102 {
103 return fluxData_;
104 }
106 const FluxData& fluxData() const
107 {
108 return fluxData_;
109 }
110
117 Scalar pressure(int phaseIdx)
118 {
119 return fluidState_->pressure(phaseIdx);
120 }
121
126 const Scalar pressure(int phaseIdx) const
127 {
128 return fluidState_->pressure(phaseIdx);
129 }
130
139 const Scalar totalConcentration(int compIdx) const
140 {
141 return massConcentration_[compIdx];
142 }
143
145 const Scalar massConcentration(int compIdx) const
146 {
147 return massConcentration_[compIdx];
148 }
149
156 void setTotalConcentration(int compIdx, Scalar value)
157 {
158 massConcentration_[compIdx] = value;
159 }
160
162 void setMassConcentration(int compIdx, Scalar value)
163 {
164 massConcentration_[compIdx] = value;
165 }
172 {
173 massConcentration_[wCompIdx] =
174 porosity * (massFraction(wPhaseIdx,wCompIdx) * saturation(wPhaseIdx) * density(wPhaseIdx)
175 + massFraction(nPhaseIdx,wCompIdx) * saturation(nPhaseIdx) * density(nPhaseIdx));
176 massConcentration_[nCompIdx] =
177 porosity * (massFraction(wPhaseIdx,nCompIdx) * saturation(wPhaseIdx) * density(wPhaseIdx)
178 + massFraction(nPhaseIdx,nCompIdx) * saturation(nPhaseIdx) * density(nPhaseIdx));
179 }
181
182
189 const Scalar& mobility(int phaseIdx) const
190 {
191 return mobility_[phaseIdx];
192 }
193
199 void setMobility(int phaseIdx, Scalar value)
200 {
201 mobility_[phaseIdx]=value;
202 }
203
210 Scalar& volumeError()
211 {
212 return volumeError_;
213 }
214
216 const Scalar& volumeError() const
217 {
218 return volumeError_;
219 }
220
228 {
229 return errorCorrection_;
230 }
231
233 const Scalar& errorCorrection() const
234 {
235 return errorCorrection_;
236 }
237
242 Scalar& dv_dp()
243 {
244 return dv_dp_;
245 }
246
248 const Scalar& dv_dp() const
249 {
250 return dv_dp_;
251 }
252
258 Scalar& dv(int compIdx)
259 {
260 return dv_[compIdx];
261 }
262
264 const Scalar& dv(int compIdx) const
265 {
266 return dv_[compIdx];
267 }
268
274 Scalar& perimeter()
275 {
276 return perimeter_;
277 }
279 const Scalar& perimeter() const
280 {
281 return perimeter_;
282 }
283
284 /*** b) from fluidstate ***/
285
287 const Scalar saturation(int phaseIdx) const
288 {
289 return fluidState_->saturation(phaseIdx);
290 }
291
293 const Scalar viscosity(int phaseIdx) const
294 {
295 return fluidState_->viscosity(phaseIdx);
296 }
297
301 const Scalar capillaryPressure() const
302 {
303 return fluidState_->pressure(nPhaseIdx) - fluidState_->pressure(wPhaseIdx);
304 }
305
307 const Scalar density(int phaseIdx) const
308 {
309 return (fluidState_->density(phaseIdx));
310 }
311
313 const Scalar massFraction(int phaseIdx, int compIdx) const
314 {
315 return fluidState_->massFraction(phaseIdx, compIdx);
316 }
317
319 const Scalar moleFraction(int phaseIdx, int compIdx) const
320 {
321 return fluidState_->moleFraction(phaseIdx, compIdx);
322 }
323
325 const Scalar temperature(int phaseIdx) const
326 {
327 return fluidState_->temperature(phaseIdx);
328 }
329
331 const Scalar phaseMassFraction(int phaseIdx) const
332 {
333 return fluidState_->phaseMassFraction(phaseIdx);
334 }
336
338 const FluidState& fluidState() const
339 {
340 return *fluidState_;
341 }
342
349 {
350 if(!fluidState_)
351 fluidState_ = std::make_shared<FluidState>();
352 return *fluidState_;
353 }
354
357 { return globalIdx_;}
368 void reset()
369 {
370 wasRefined_ = false;
372 }
375 { return wasRefined_;}
377 const bool& wasRefined() const
378 { return wasRefined_;}
379
385 const bool& isUpwindCell(int indexInInside, int phaseIdx) const
386 {
387 return fluxData_.isUpwindCell(indexInInside, phaseIdx);
388 }
395 void setUpwindCell(int indexInInside, int phaseIdx, bool value)
396 {
397 fluxData_.setUpwindCell(indexInInside, phaseIdx, value);
398 }
399
400};
401}
402#endif
Definition: adapt.hh:29
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:149
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
Storage container for discretized data of the constitutive relations for one element.
Definition: 2p2c/sequential/celldata.hh:44
Scalar volumeError_
Definition: 2p2c/sequential/celldata.hh:70
int globalIdx_
Definition: 2p2c/sequential/celldata.hh:77
void volumeDerivativesAvailable(bool value)
Specifies if volume derivatives are computed and available.
Definition: 2p2c/sequential/celldata.hh:365
FluidState & manipulateFluidState()
Allows manipulation of the cells fluid state Fluidstate is stored as a pointer, initialized as a null...
Definition: 2p2c/sequential/celldata.hh:348
FluxData fluxData_
Definition: 2p2c/sequential/celldata.hh:81
void setTotalConcentration(int compIdx, Scalar value)
Sets the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:156
bool hasVolumeDerivatives() const
Indicates if volume derivatives are computed and available.
Definition: 2p2c/sequential/celldata.hh:359
void setUpwindCell(int indexInInside, int phaseIdx, bool value)
Specifies if current cell is the upwind cell for a given interface.
Definition: 2p2c/sequential/celldata.hh:395
Scalar mobility_[numPhases]
Definition: 2p2c/sequential/celldata.hh:67
const Scalar viscosity(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:293
Scalar & volumeError()
Return the volume error [-]. This quantity stands for the deviation of real fluid volume to available...
Definition: 2p2c/sequential/celldata.hh:210
void confirmVolumeDerivatives()
Specifies that volume derivatives are computed and available.
Definition: 2p2c/sequential/celldata.hh:362
const Scalar temperature(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:325
const Scalar & volumeError() const
Return the volume error [-].
Definition: 2p2c/sequential/celldata.hh:216
bool & wasRefined()
Indicates if current cell was refined at this time step.
Definition: 2p2c/sequential/celldata.hh:374
CellData2P2C()
Constructor for a local CellData object.
Definition: 2p2c/sequential/celldata.hh:85
Scalar dv_[numComponents]
Definition: 2p2c/sequential/celldata.hh:73
bool volumeDerivativesAvailable_
Definition: 2p2c/sequential/celldata.hh:74
const Scalar saturation(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:287
const Scalar massConcentration(int compIdx) const
Returns the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:145
int & globalIdx()
stores this cell datas index, only for debugging purposes!!
Definition: 2p2c/sequential/celldata.hh:356
Scalar perimeter_
Definition: 2p2c/sequential/celldata.hh:78
void setMassConcentration(int compIdx, Scalar value)
Sets the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:162
Scalar pressure(int phaseIdx)
Acess to the phase pressure.
Definition: 2p2c/sequential/celldata.hh:117
Scalar dv_dp_
Definition: 2p2c/sequential/celldata.hh:72
void reset()
Resets the cell data after a timestep was completed: No volume derivatives yet available.
Definition: 2p2c/sequential/celldata.hh:368
std::shared_ptr< FluidState > fluidState_
Definition: 2p2c/sequential/celldata.hh:80
const Scalar density(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:307
const Scalar totalConcentration(int compIdx) const
Returns the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:139
Scalar errorCorrection_
Definition: 2p2c/sequential/celldata.hh:71
void setMobility(int phaseIdx, Scalar value)
Set phase mobilities.
Definition: 2p2c/sequential/celldata.hh:199
const Scalar & dv(int compIdx) const
Return the derivative of spec. volume w.r.t. change of mass For details, see description of FVPressur...
Definition: 2p2c/sequential/celldata.hh:264
const FluidState & fluidState() const
Returns a reference to the cells fluid state.
Definition: 2p2c/sequential/celldata.hh:338
void calculateMassConcentration(Scalar porosity)
Calculate the total mass concentration of a component for a given porosity (within the initializatio...
Definition: 2p2c/sequential/celldata.hh:171
const Scalar & errorCorrection() const
Return the error Correction.
Definition: 2p2c/sequential/celldata.hh:233
const Scalar phaseMassFraction(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:331
Scalar & dv_dp()
Return the derivative of specific volume w.r.t. pressure For details, see description of FVPressureCo...
Definition: 2p2c/sequential/celldata.hh:242
Scalar & errorCorrection()
Return the error Correction This quantifies the damped error that actually entered the pressure equat...
Definition: 2p2c/sequential/celldata.hh:227
const bool & isUpwindCell(int indexInInside, int phaseIdx) const
Indicates if current cell is the upwind cell for a given interface.
Definition: 2p2c/sequential/celldata.hh:385
FluxData & fluxData()
Acess to flux data, representing information living on the intersections.
Definition: 2p2c/sequential/celldata.hh:101
Scalar & dv(int compIdx)
Return the derivative of spec. volume w.r.t. change of mass For details, see description of FVPressur...
Definition: 2p2c/sequential/celldata.hh:258
const Scalar moleFraction(int phaseIdx, int compIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:319
const Scalar & mobility(int phaseIdx) const
Return phase mobilities.
Definition: 2p2c/sequential/celldata.hh:189
const Scalar & dv_dp() const
Return the derivative of specific volume w.r.t. pressure.
Definition: 2p2c/sequential/celldata.hh:248
Scalar massConcentration_[numComponents]
Definition: 2p2c/sequential/celldata.hh:65
Scalar & perimeter()
Return cell perimeter (as weithing function) The cell perimeter is used in combination with the face ...
Definition: 2p2c/sequential/celldata.hh:274
const Scalar massFraction(int phaseIdx, int compIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:313
const Scalar & perimeter() const
Return cell perimeter (as weithing function)
Definition: 2p2c/sequential/celldata.hh:279
const Scalar pressure(int phaseIdx) const
Acess to the phase pressure.
Definition: 2p2c/sequential/celldata.hh:126
const FluxData & fluxData() const
Constant acess to flux data, representing information living on the intersections.
Definition: 2p2c/sequential/celldata.hh:106
bool wasRefined_
Definition: 2p2c/sequential/celldata.hh:75
const bool & wasRefined() const
Indicates if current cell was refined at this time step.
Definition: 2p2c/sequential/celldata.hh:377
const Scalar capillaryPressure() const
Returns the capillary pressure .
Definition: 2p2c/sequential/celldata.hh:301
Class including the variables and data of discretized data of the constitutive relations.
Definition: 2p2c/sequential/fluxdata.hh:43
void setUpwindCell(int indexInInside, int equationIdx, bool value)
Sets the upwind information.
Definition: 2p2c/sequential/fluxdata.hh:170
const bool & isUpwindCell(int indexInInside, int equationIdx) const
Functions returning upwind information.
Definition: 2p2c/sequential/fluxdata.hh:159
Class including the variables and data of discretized data of the constitutive relations.