3.1-git
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:
46 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
48 using FluidState = typename GET_PROP_TYPE(TypeTag, FluidState);
49
50 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
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 = GET_PROP_VALUE(TypeTag, NumPhases),
61 numComponents = GET_PROP_VALUE(TypeTag, 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
136 [[deprecated("cellData.setPressure is deprecated and will be removed after release 3.1. Use cellData.manipulateFluidState() and set pressure in fluid state directly")]]
137 void setPressure(int phaseIdx, Scalar value)
138 {
139 manipulateFluidState().setPressure(phaseIdx, value);
140 }
141
150 const Scalar totalConcentration(int compIdx) const
151 {
152 return massConcentration_[compIdx];
153 }
154
156 const Scalar massConcentration(int compIdx) const
157 {
158 return massConcentration_[compIdx];
159 }
160
167 void setTotalConcentration(int compIdx, Scalar value)
168 {
169 massConcentration_[compIdx] = value;
170 }
171
173 void setMassConcentration(int compIdx, Scalar value)
174 {
175 massConcentration_[compIdx] = value;
176 }
183 {
184 massConcentration_[wCompIdx] =
185 porosity * (massFraction(wPhaseIdx,wCompIdx) * saturation(wPhaseIdx) * density(wPhaseIdx)
186 + massFraction(nPhaseIdx,wCompIdx) * saturation(nPhaseIdx) * density(nPhaseIdx));
187 massConcentration_[nCompIdx] =
188 porosity * (massFraction(wPhaseIdx,nCompIdx) * saturation(wPhaseIdx) * density(wPhaseIdx)
189 + massFraction(nPhaseIdx,nCompIdx) * saturation(nPhaseIdx) * density(nPhaseIdx));
190 }
192
193
200 const Scalar& mobility(int phaseIdx) const
201 {
202 return mobility_[phaseIdx];
203 }
204
210 void setMobility(int phaseIdx, Scalar value)
211 {
212 mobility_[phaseIdx]=value;
213 }
214
221 Scalar& volumeError()
222 {
223 return volumeError_;
224 }
225
227 const Scalar& volumeError() const
228 {
229 return volumeError_;
230 }
231
239 {
240 return errorCorrection_;
241 }
242
244 const Scalar& errorCorrection() const
245 {
246 return errorCorrection_;
247 }
248
253 Scalar& dv_dp()
254 {
255 return dv_dp_;
256 }
257
259 const Scalar& dv_dp() const
260 {
261 return dv_dp_;
262 }
263
269 Scalar& dv(int compIdx)
270 {
271 return dv_[compIdx];
272 }
273
275 const Scalar& dv(int compIdx) const
276 {
277 return dv_[compIdx];
278 }
279
285 Scalar& perimeter()
286 {
287 return perimeter_;
288 }
290 const Scalar& perimeter() const
291 {
292 return perimeter_;
293 }
294
295 /*** b) from fluidstate ***/
296
298 [[deprecated("cellData.setSaturation is deprecated and will be removed after release 3.1. Use cellData.manipulateFluidState() and set saturation in fluid state directly")]]
299 void setSaturation(int phaseIdx, Scalar value)
300 {
301 fluidState_->setSaturation(phaseIdx, value);
302 fluidState_->setSaturation(1-phaseIdx, 1.0-value);
303 }
304
306 const Scalar saturation(int phaseIdx) const
307 {
308 return fluidState_->saturation(phaseIdx);
309 }
310
312 [[deprecated("cellData.setViscosity is deprecated and will be removed after release 3.1. Use cellData.manipulateFluidState() and set viscosity in fluid state directly")]]
313 void setViscosity(int phaseIdx, Scalar value)
314 {
315 fluidState_->setViscosity(phaseIdx, value);
316 }
318 const Scalar viscosity(int phaseIdx) const
319 {
320 return fluidState_->viscosity(phaseIdx);
321 }
322
326 const Scalar capillaryPressure() const
327 {
328 return fluidState_->pressure(nPhaseIdx) - fluidState_->pressure(wPhaseIdx);
329 }
330
332 const Scalar density(int phaseIdx) const
333 {
334 return (fluidState_->density(phaseIdx));
335 }
336
338 const Scalar massFraction(int phaseIdx, int compIdx) const
339 {
340 return fluidState_->massFraction(phaseIdx, compIdx);
341 }
342
344 const Scalar moleFraction(int phaseIdx, int compIdx) const
345 {
346 return fluidState_->moleFraction(phaseIdx, compIdx);
347 }
348
350 const Scalar temperature(int phaseIdx) const
351 {
352 return fluidState_->temperature(phaseIdx);
353 }
354
356 const Scalar phaseMassFraction(int phaseIdx) const
357 {
358 return fluidState_->phaseMassFraction(phaseIdx);
359 }
361
363 const FluidState& fluidState() const
364 {
365 return *fluidState_;
366 }
367
374 {
375 if(!fluidState_)
376 fluidState_ = std::make_shared<FluidState>();
377 return *fluidState_;
378 }
379
382 { return globalIdx_;}
393 void reset()
394 {
395 wasRefined_ = false;
397 }
400 { return wasRefined_;}
402 const bool& wasRefined() const
403 { return wasRefined_;}
404
410 const bool& isUpwindCell(int indexInInside, int phaseIdx) const
411 {
412 return fluxData_.isUpwindCell(indexInInside, phaseIdx);
413 }
420 void setUpwindCell(int indexInInside, int phaseIdx, bool value)
421 {
422 fluxData_.setUpwindCell(indexInInside, phaseIdx, value);
423 }
424
425};
426}
427#endif
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag NumComponents
Number of components in the system.
Definition: porousmediumflow/sequential/properties.hh:70
Property tag NumPhases
Number of phases in the system.
Definition: porousmediumflow/sequential/properties.hh:69
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:390
FluidState & manipulateFluidState()
Allows manipulation of the cells fluid state Fluidstate is stored as a pointer, initialized as a null...
Definition: 2p2c/sequential/celldata.hh:373
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:167
bool hasVolumeDerivatives() const
Indicates if volume derivatives are computed and available.
Definition: 2p2c/sequential/celldata.hh:384
void setViscosity(int phaseIdx, Scalar value)
DOC ME!
Definition: 2p2c/sequential/celldata.hh:313
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:420
Scalar mobility_[numPhases]
Definition: 2p2c/sequential/celldata.hh:67
void setPressure(int phaseIdx, Scalar value)
Modify the phase pressure.
Definition: 2p2c/sequential/celldata.hh:137
const Scalar viscosity(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:318
Scalar & volumeError()
Return the volume error [-]. This quantity stands for the deviation of real fluid volume to available...
Definition: 2p2c/sequential/celldata.hh:221
void confirmVolumeDerivatives()
Specifies that volume derivatives are computed and available.
Definition: 2p2c/sequential/celldata.hh:387
const Scalar temperature(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:350
const Scalar & volumeError() const
Return the volume error [-].
Definition: 2p2c/sequential/celldata.hh:227
bool & wasRefined()
Indicates if current cell was refined at this time step.
Definition: 2p2c/sequential/celldata.hh:399
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:306
const Scalar massConcentration(int compIdx) const
Returns the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:156
int & globalIdx()
stores this cell datas index, only for debugging purposes!!
Definition: 2p2c/sequential/celldata.hh:381
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:173
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:393
std::shared_ptr< FluidState > fluidState_
Definition: 2p2c/sequential/celldata.hh:80
const Scalar density(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:332
const Scalar totalConcentration(int compIdx) const
Returns the total mass concentration of a component .
Definition: 2p2c/sequential/celldata.hh:150
Scalar errorCorrection_
Definition: 2p2c/sequential/celldata.hh:71
void setMobility(int phaseIdx, Scalar value)
Set phase mobilities.
Definition: 2p2c/sequential/celldata.hh:210
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:275
const FluidState & fluidState() const
Returns a reference to the cells fluid state.
Definition: 2p2c/sequential/celldata.hh:363
void calculateMassConcentration(Scalar porosity)
Calculate the total mass concentration of a component for a given porosity (within the initializatio...
Definition: 2p2c/sequential/celldata.hh:182
const Scalar & errorCorrection() const
Return the error Correction.
Definition: 2p2c/sequential/celldata.hh:244
const Scalar phaseMassFraction(int phaseIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:356
Scalar & dv_dp()
Return the derivative of specific volume w.r.t. pressure For details, see description of FVPressureCo...
Definition: 2p2c/sequential/celldata.hh:253
Scalar & errorCorrection()
Return the error Correction This quantifies the damped error that actually entered the pressure equat...
Definition: 2p2c/sequential/celldata.hh:238
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:410
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:269
const Scalar moleFraction(int phaseIdx, int compIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:344
const Scalar & mobility(int phaseIdx) const
Return phase mobilities.
Definition: 2p2c/sequential/celldata.hh:200
const Scalar & dv_dp() const
Return the derivative of specific volume w.r.t. pressure.
Definition: 2p2c/sequential/celldata.hh:259
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:285
const Scalar massFraction(int phaseIdx, int compIdx) const
DOC ME!
Definition: 2p2c/sequential/celldata.hh:338
const Scalar & perimeter() const
Return cell perimeter (as weithing function)
Definition: 2p2c/sequential/celldata.hh:290
void setSaturation(int phaseIdx, Scalar value)
DOC ME!
Definition: 2p2c/sequential/celldata.hh:299
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:402
const Scalar capillaryPressure() const
Returns the capillary pressure .
Definition: 2p2c/sequential/celldata.hh:326
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.