My Project
Loading...
Searching...
No Matches
Aquancon.hpp
1/*
2 Copyright (C) 2017 TNO
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 3 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
20#ifndef OPM_AQUANCON_HPP
21#define OPM_AQUANCON_HPP
22
23/*
24 Aquancon is a data container object meant to hold the data from the AQUANCON keyword.
25 This also includes the logic for parsing and connections to grid cells. It is meant to be used by opm-grid and opm-simulators in order to
26 implement the analytical aquifer models in OPM Flow.
27*/
28
29#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
30
31#include <cstddef>
32#include <unordered_map>
33#include <vector>
34
35namespace Opm {
36 class EclipseGrid;
37 class Deck;
38}
39
40namespace Opm { namespace RestartIO {
41 class RstAquifer;
42}} // Opm::RestartIO
43
44namespace Opm {
45
46 class Aquancon {
47 public:
48
49 struct AquancCell {
50 int aquiferID;
51 std::size_t global_index;
52 double influx_coeff;
53 double effective_facearea; // Needed for restart output only.
54 FaceDir::DirEnum face_dir;
55
56 AquancCell(const int aquiferID_arg,
57 const std::size_t gi,
58 const double ic,
59 const double eff_faceArea,
60 const FaceDir::DirEnum fd)
61 : aquiferID(aquiferID_arg)
62 , global_index(gi)
63 , influx_coeff(ic)
64 , effective_facearea(eff_faceArea)
65 , face_dir(fd)
66 {}
67
68 AquancCell() = default;
69
70 bool operator==(const AquancCell& other) const {
71 return (this->aquiferID == other.aquiferID)
72 && (this->global_index == other.global_index)
73 && (this->influx_coeff == other.influx_coeff)
74 && (this->effective_facearea == other.effective_facearea)
75 && (this->face_dir == other.face_dir);
76 }
77
78 template<class Serializer>
79 void serializeOp(Serializer& serializer)
80 {
81 serializer(this->aquiferID);
82 serializer(this->global_index);
83 serializer(this->influx_coeff);
84 serializer(this->effective_facearea);
85 serializer(this->face_dir);
86 }
87 };
88
89 Aquancon() = default;
90 Aquancon(const EclipseGrid& grid, const Deck& deck);
91 explicit Aquancon(const std::unordered_map<int, std::vector<Aquancon::AquancCell>>& data);
92
93 void pruneDeactivatedAquiferConnections(const std::vector<std::size_t>& deactivated_cells);
94 void loadFromRestart(const RestartIO::RstAquifer& rst_aquifers);
95
96 static Aquancon serializationTestObject();
97
98 const std::unordered_map<int, std::vector<Aquancon::AquancCell>>& data() const;
99 bool operator==(const Aquancon& other) const;
100 bool active() const;
101
102 bool hasAquiferConnections(int aquiferID) const;
103 const std::vector<Aquancon::AquancCell>& getConnections(int aquiferID) const;
104
105 template<class Serializer>
106 void serializeOp(Serializer& serializer)
107 {
108 serializer(this->cells);
109 }
110
111 private:
112 std::unordered_map<int, std::vector<Aquancon::AquancCell>> cells;
113 };
114}
115
116#endif
Definition Aquancon.hpp:46
Definition Deck.hpp:49
About cell information and dimension: The actual grid information is held in a pointer to an ERT ecl_...
Definition EclipseGrid.hpp:55
Definition aquifer.hpp:45
Class for (de-)serializing.
Definition Serializer.hpp:84
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Aquancon.hpp:49