Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup Components
* \brief A generic template for tabulated material laws that depend
* on two parameters.
*/
#ifndef DUMUX_EXERCISE_CO2TABLES_HH
#define DUMUX_EXERCISE_CO2TABLES_HH
#include <dune/common/float_cmp.hh>
namespace Dumux::BioMinCO2Tables {
/*!
* \ingroup Components
* \brief A generic template for tabulated material laws that depend
* on two parameters.
*/
template <class Traits>
class TabulatedProperties
{
using Scalar = typename Traits::Scalar;
static constexpr auto numTempSteps = Traits::numTempSteps;
static constexpr auto numPressSteps = Traits::numPressSteps;
public:
TabulatedProperties() = default;
constexpr Scalar minTemp() const { return Traits::minTemp; }
constexpr Scalar maxTemp() const { return Traits::maxTemp; }
constexpr Scalar minPress() const { return Traits::minPress; }
constexpr Scalar maxPress() const { return Traits::maxPress; }
constexpr bool applies(Scalar temperature, Scalar pressure) const
{
return minTemp() <= temperature && temperature <= maxTemp() &&
minPress() <= pressure && pressure <= maxPress();
}
constexpr Scalar at(Scalar temperature, Scalar pressure) const
{
if (!applies(temperature, pressure))
{
if (temperature<minTemp()) temperature = minTemp();
else if (temperature>maxTemp()) temperature = maxTemp();
if (pressure<minPress()) pressure = minPress();
else if (pressure>maxPress()) pressure = maxPress();
}
const int i = findTempIdx_(temperature);
const int j = findPressIdx_(pressure);
const Scalar tempAtI = temperatureAt_(i);
const Scalar tempAtI1 = temperatureAt_(i + 1);
const Scalar pressAtI = pressureAt_(j);
const Scalar pressAtI1 = pressureAt_(j + 1);
const Scalar alpha = (temperature - tempAtI)/(tempAtI1 - tempAtI);
const Scalar beta = (pressure - pressAtI)/(pressAtI1 - pressAtI);
// bi-linear interpolation
const Scalar lowresValue =
(1-alpha)*(1-beta)*val(i, j) +
(1-alpha)*( beta)*val(i, j + 1) +
( alpha)*(1-beta)*val(i + 1, j) +
( alpha)*( beta)*val(i + 1, j + 1);
// return the weighted sum of the low- and high-resolution values
return lowresValue;
}
constexpr Scalar val(int i, int j) const
{ return Traits::vals[i][j]; }
private:
constexpr int findTempIdx_(Scalar temperature) const
{
if (Dune::FloatCmp::eq<Scalar>(temperature, maxTemp()))
return numTempSteps - 2;
const int result = static_cast<int>((temperature - minTemp())/(maxTemp() - minTemp())*(numTempSteps - 1));
using std::clamp;
return clamp(result, 0, numTempSteps - 2);
}
constexpr int findPressIdx_(Scalar pressure) const
{
if (Dune::FloatCmp::eq<Scalar>(pressure, maxPress()))
return numPressSteps - 2;
const int result = static_cast<int>((pressure - minPress())/(maxPress() - minPress())*(numPressSteps - 1));
using std::clamp;
return clamp(result, 0, numPressSteps - 2);
}
constexpr Scalar temperatureAt_(int i) const
{ return i*(maxTemp() - minTemp())/(numTempSteps - 1) + minTemp(); }
constexpr Scalar pressureAt_(int j) const
{ return j*(maxPress() - minPress())/(numPressSteps - 1) + minPress(); }
};
#ifndef DOXYGEN // hide from doxygen
// the real work is done by some external program which provides
// ready-to-use tables.
#include "co2values.inc"
#endif
using TabulatedDensity = TabulatedProperties<TabulatedDensityTraits>;
using TabulatedEnthalpy = TabulatedProperties<TabulatedEnthalpyTraits>;
// this class collects all the tabulated quantities in one convenient place
struct CO2Tables
{
static constexpr inline TabulatedEnthalpy tabulatedEnthalpy = {};
static constexpr inline TabulatedDensity tabulatedDensity = {};
};
} // end namespace Dumux::GeneratedCO2Tables
#endif