20#ifndef OPM_OUTPUT_WELLS_HPP
21#define OPM_OUTPUT_WELLS_HPP
23#include <opm/common/ErrorMacros.hpp>
24#include <opm/common/OpmLog/OpmLog.hpp>
25#include <opm/output/data/GuideRateValue.hpp>
26#include <opm/input/eclipse/Schedule/Well/WellEnums.hpp>
36#include <unordered_map>
39namespace Opm {
namespace data {
53 enum class opt : uint32_t {
60 dissolved_gas = (1 << 6),
61 vaporized_oil = (1 << 7),
62 reservoir_water = (1 << 8),
63 reservoir_oil = (1 << 9),
64 reservoir_gas = (1 << 10),
65 productivity_index_water = (1 << 11),
66 productivity_index_oil = (1 << 12),
67 productivity_index_gas = (1 << 13),
68 well_potential_water = (1 << 14),
69 well_potential_oil = (1 << 15),
70 well_potential_gas = (1 << 16),
74 microbial = (1 << 20),
77 vaporized_water = (1 << 23),
83 using enum_size = std::underlying_type< opt >::type;
86 inline bool has( opt )
const;
90 inline double get( opt m )
const;
93 inline double get( opt m,
double default_value )
const;
94 inline double get( opt m,
double default_value ,
const std::string& tracer_name )
const;
98 inline Rates&
set( opt m,
double value );
99 inline Rates&
set( opt m,
double value ,
const std::string& tracer_name );
104 template <
class MessageBufferType>
105 void write(MessageBufferType& buffer)
const;
106 template <
class MessageBufferType>
107 void read(MessageBufferType& buffer);
109 bool operator==(
const Rates& rat2)
const;
111 template<
class Serializer>
121 serializer(dissolved_gas);
122 serializer(vaporized_oil);
123 serializer(reservoir_water);
124 serializer(reservoir_oil);
125 serializer(reservoir_gas);
126 serializer(productivity_index_water);
127 serializer(productivity_index_oil);
128 serializer(productivity_index_gas);
129 serializer(well_potential_water);
130 serializer(well_potential_oil);
131 serializer(well_potential_gas);
135 serializer(microbial);
138 serializer(vaporized_water);
139 serializer(mass_gas);
140 serializer(mass_wat);
141 serializer(wat_frac);
144 static Rates serializationTestObject()
147 rat1.
set(opt::wat, 1.0);
148 rat1.
set(opt::oil, 2.0);
149 rat1.
set(opt::gas, 3.0);
150 rat1.
set(opt::polymer, 4.0);
151 rat1.
set(opt::solvent, 5.0);
152 rat1.
set(opt::energy, 6.0);
153 rat1.
set(opt::dissolved_gas, 7.0);
154 rat1.
set(opt::vaporized_oil, 8.0);
155 rat1.
set(opt::reservoir_water, 9.0);
156 rat1.
set(opt::reservoir_oil, 10.0);
157 rat1.
set(opt::reservoir_gas, 11.0);
158 rat1.
set(opt::productivity_index_water, 12.0);
159 rat1.
set(opt::productivity_index_oil, 13.0);
160 rat1.
set(opt::productivity_index_gas, 14.0);
161 rat1.
set(opt::well_potential_water, 15.0);
162 rat1.
set(opt::well_potential_oil, 16.0);
163 rat1.
set(opt::well_potential_gas, 17.0);
164 rat1.
set(opt::brine, 18.0);
165 rat1.
set(opt::alq, 19.0);
166 rat1.
set(opt::microbial, 20.0);
167 rat1.
set(opt::oxygen, 21.0);
168 rat1.
set(opt::urea, 22.0);
169 rat1.
set(opt::vaporized_water, 23.0);
170 rat1.
set(opt::mass_gas, 24.0);
171 rat1.
set(opt::mass_wat, 25.0);
172 rat1.
set(opt::wat_frac, 26.0);
173 rat1.tracer.insert({
"test_tracer", 1.0});
179 double& get_ref( opt );
180 double& get_ref( opt,
const std::string& tracer_name );
181 const double& get_ref( opt )
const;
182 const double& get_ref( opt,
const std::string& tracer_name )
const;
184 opt mask =
static_cast< opt
>( 0 );
189 double polymer = 0.0;
190 double solvent = 0.0;
192 double dissolved_gas = 0.0;
193 double vaporized_oil = 0.0;
194 double reservoir_water = 0.0;
195 double reservoir_oil = 0.0;
196 double reservoir_gas = 0.0;
197 double productivity_index_water = 0.0;
198 double productivity_index_oil = 0.0;
199 double productivity_index_gas = 0.0;
200 double well_potential_water = 0.0;
201 double well_potential_oil = 0.0;
202 double well_potential_gas = 0.0;
205 std::map<std::string, double> tracer{};
206 double microbial = 0.0;
209 double vaporized_water = 0.0;
210 double mass_gas = 0.0;
211 double mass_wat = 0.0;
212 double wat_frac = 0.0;
219 double skin_factor{};
224 double area_of_flow{};
225 double flow_factor{};
226 double fracture_rate{};
228 template <
class Serializer>
233 serializer(skin_factor);
234 serializer(thickness);
238 serializer(area_of_flow);
239 serializer(flow_factor);
240 serializer(fracture_rate);
245 return (this->rate == filtrate.rate)
246 && (this->total == filtrate.total)
247 && (this->skin_factor == filtrate.skin_factor)
248 && (this->thickness == filtrate.thickness)
249 && (this->perm == filtrate.perm)
250 && (this->poro == filtrate.poro)
251 && (this->radius == filtrate.radius)
252 && (this->area_of_flow == filtrate.area_of_flow)
253 && (this->flow_factor == filtrate.flow_factor)
254 && (this->fracture_rate == filtrate.fracture_rate)
260 return {0.8, 100.0, -1.0, 2.0, 1.0e-9,
261 0.3, 0.05, 0.8, 0.1, 0.7};
264 template <
class MessageBufferType>
265 void write(MessageBufferType& buffer)
const;
267 template <
class MessageBufferType>
268 void read(MessageBufferType& buffer);
279 double filter_volume{};
281 double avg_filter_width{};
282 double inj_pressure{};
284 double inj_wellrate{};
286 template <
class Serializer>
295 serializer(filter_volume);
296 serializer(avg_width);
297 serializer(avg_filter_width);
298 serializer(inj_pressure);
300 serializer(inj_wellrate);
305 return (this->area == fraccon.area)
306 && (this->flux == fraccon.flux)
307 && (this->height == fraccon.height)
308 && (this->length == fraccon.length)
309 && (this->WI == fraccon.WI)
310 && (this->volume == fraccon.volume)
311 && (this->filter_volume == fraccon.filter_volume)
312 && (this->avg_width == fraccon.avg_width)
313 && (this->avg_filter_width == fraccon.avg_filter_width)
314 && (this->inj_pressure == fraccon.inj_pressure)
315 && (this->inj_bhp == fraccon.inj_bhp)
316 && (this->inj_wellrate == fraccon.inj_wellrate)
322 return {0.8, 100.0, 1.3, 1.4, 10.0, 4.0, 0.4, 0.5, 0.05, 200.0, 150.0, 500.0};
325 template <
class MessageBufferType>
326 void write(MessageBufferType& buffer)
const
328 buffer.write(this->area);
329 buffer.write(this->flux);
330 buffer.write(this->height);
331 buffer.write(this->length);
332 buffer.write(this->WI);
333 buffer.write(this->volume);
334 buffer.write(this->filter_volume);
335 buffer.write(this->avg_width);
336 buffer.write(this->avg_filter_width);
337 buffer.write(this->inj_pressure);
338 buffer.write(this->inj_bhp);
339 buffer.write(this->inj_wellrate);
342 template <
class MessageBufferType>
343 void read(MessageBufferType& buffer)
345 buffer.read(this->area);
346 buffer.read(this->flux);
347 buffer.read(this->height);
348 buffer.read(this->length);
349 buffer.read(this->WI);
350 buffer.read(this->volume);
351 buffer.read(this->filter_volume);
352 buffer.read(this->avg_width);
353 buffer.read(this->avg_filter_width);
354 buffer.read(this->inj_pressure);
355 buffer.read(this->inj_bhp);
356 buffer.read(this->inj_wellrate);
384 12.34, 56.78, 9.10, 11.12
393 template <
class Serializer>
396 serializer(this->avg);
397 serializer(this->max);
398 serializer(this->min);
399 serializer(this->stdev);
411 return (this->avg == that.
avg)
412 && (this->max == that.
max)
413 && (this->min == that.
min)
414 && (this->stdev == that.
stdev)
419 template <
class MessageBufferType>
420 void write(MessageBufferType& buffer)
const
422 buffer.write(this->avg);
423 buffer.write(this->max);
424 buffer.write(this->min);
425 buffer.write(this->stdev);
429 template <
class MessageBufferType>
430 void read(MessageBufferType& buffer)
432 buffer.read(this->avg);
433 buffer.read(this->max);
434 buffer.read(this->min);
435 buffer.read(this->stdev);
471 template <
class Serializer>
474 serializer(this->numCells);
475 serializer(this->press);
476 serializer(this->rate);
477 serializer(this->width);
489 return (this->numCells == that.
numCells)
490 && (this->press == that.
press)
491 && (this->rate == that.
rate)
492 && (this->width == that.
width)
497 template <
class MessageBufferType>
498 void write(MessageBufferType& buffer)
const
500 buffer.write(this->numCells);
501 buffer.write(this->press);
502 buffer.write(this->rate);
503 buffer.write(this->width);
507 template <
class MessageBufferType>
508 void read(MessageBufferType& buffer)
510 buffer.read(this->numCells);
511 buffer.read(this->press);
512 buffer.read(this->rate);
513 buffer.read(this->width);
519 using global_index = std::size_t;
520 static const constexpr int restart_size = 6;
522 global_index index{};
525 double reservoir_rate{};
526 double cell_pressure{};
527 double cell_saturation_water{};
528 double cell_saturation_gas{};
529 double effective_Kh{};
530 double trans_factor{};
532 double compact_mult{1.0};
543 bool operator==(
const Connection& conn2)
const
545 return (index == conn2.index)
546 && (rates == conn2.rates)
547 && (pressure == conn2.pressure)
548 && (reservoir_rate == conn2.reservoir_rate)
549 && (cell_pressure == conn2.cell_pressure)
550 && (cell_saturation_water == conn2.cell_saturation_water)
551 && (cell_saturation_gas == conn2.cell_saturation_gas)
552 && (effective_Kh == conn2.effective_Kh)
553 && (trans_factor == conn2.trans_factor)
554 && (d_factor == conn2.d_factor)
555 && (compact_mult == conn2.compact_mult)
556 && (lgr_grid == conn2.lgr_grid)
557 && (filtrate == conn2.filtrate)
558 && (fracture == conn2.fracture)
559 && (this->fract == conn2.
fract)
563 template <
class MessageBufferType>
564 void write(MessageBufferType& buffer)
const;
565 template <
class MessageBufferType>
566 void read(MessageBufferType& buffer);
568 template<
class Serializer>
573 serializer(pressure);
574 serializer(reservoir_rate);
575 serializer(cell_pressure);
576 serializer(cell_saturation_water);
577 serializer(cell_saturation_gas);
578 serializer(effective_Kh);
579 serializer(trans_factor);
580 serializer(d_factor);
581 serializer(compact_mult);
582 serializer(lgr_grid);
583 serializer(filtrate);
584 serializer(fracture);
585 serializer(this->fract);
588 static Connection serializationTestObject()
591 1, Rates::serializationTestObject(),
593 6.0, 7.0, 8.0, 9.0, 0.987,
595 ConnectionFiltrate::serializationTestObject(),
596 ConnectionFracture::serializationTestObject(),
605 enum class Value : std::size_t {
606 Pressure, PDrop, PDropHydrostatic, PDropAccel, PDropFriction,
609 double& operator[](
const Value i)
611 return this->values_[this->index(i)];
614 double operator[](
const Value i)
const
616 return this->values_[this->index(i)];
621 return this->values_ == segpres2.values_;
624 template <
class MessageBufferType>
625 void write(MessageBufferType& buffer)
const
627 for (
const auto& value : this->values_) {
632 template <
class MessageBufferType>
633 void read(MessageBufferType& buffer)
635 for (
auto& value : this->values_) {
640 template<
class Serializer>
649 spres[Value::Pressure] = 1.0;
650 spres[Value::PDrop] = 2.0;
651 spres[Value::PDropHydrostatic] = 3.0;
652 spres[Value::PDropAccel] = 4.0;
653 spres[Value::PDropFriction] = 5.0;
659 constexpr static std::size_t numvals = 5;
661 std::array<double, numvals> values_ = {0};
663 std::size_t index(
const Value ix)
const
665 return static_cast<std::size_t
>(ix);
669 template <
typename Items>
673 using Item =
typename Items::Item;
677 this->has_ =
static_cast<unsigned char>(0);
678 this->value_.fill(0.0);
681 constexpr bool has(
const Item p)
const
683 const auto i = this->index(p);
685 return (i < Size) && this->hasItem(i);
690 return (this->has_ == that.has_)
691 && (this->value_ == that.value_);
694 double get(
const Item p)
const
696 if (! this->has(p)) {
697 throw std::invalid_argument {
698 "Request for Unset Item Value for " + Items::itemName(p)
702 return this->value_[ this->index(p) ];
707 const auto i = this->index(p);
710 throw std::invalid_argument {
711 "Cannot Assign Item Value for Unsupported Item '"
712 + Items::itemName(p) +
'\''
716 this->has_ |= 1 << i;
717 this->value_[i] = value;
722 template <
class MessageBufferType>
723 void write(MessageBufferType& buffer)
const
725 buffer.write(this->has_);
727 for (
const auto& x : this->value_) {
732 template <
class MessageBufferType>
733 void read(MessageBufferType& buffer)
736 buffer.read(this->has_);
738 for (
auto& x : this->value_) {
743 template <
class Serializer>
746 serializer(this->has_);
747 serializer(this->value_);
754 for (
const auto& [item, value] : Items::serializationTestItems()) {
755 quant.set(item, value);
762 enum { Size =
static_cast<std::size_t
>(Item::NumItems) };
764 static_assert(Size <= static_cast<std::size_t>(CHAR_BIT),
765 "Number of items must not exceed CHAR_BIT");
769 unsigned char has_{};
772 std::array<double, Size> value_{};
774 constexpr std::size_t index(
const Item p)
const noexcept
776 return static_cast<std::size_t
>(p);
779 bool hasItem(
const std::size_t i)
const
781 return (this->has_ & (1 << i)) != 0;
794 static std::string itemName(
const Item p)
797 case Item::Oil:
return "Oil";
798 case Item::Gas:
return "Gas";
799 case Item::Water:
return "Water";
802 return "Out of bounds (NumItems)";
805 return "Unknown (" + std::to_string(
static_cast<int>(p)) +
')';
808 static auto serializationTestItems()
811 std::pair { Item::Oil , 1.0 },
812 std::pair { Item::Gas , 7.0 },
813 std::pair { Item::Water, 2.9 },
821 Oil, Gas, Water, Mixture, MixtureWithExponents,
827 static std::string itemName(
const Item p)
830 case Item::Oil:
return "Oil";
831 case Item::Gas:
return "Gas";
832 case Item::Water:
return "Water";
833 case Item::Mixture:
return "Mixture";
834 case Item::MixtureWithExponents:
return "MixtureWithExponents";
837 return "Out of bounds (NumItems)";
840 return "Unknown (" + std::to_string(
static_cast<int>(p)) +
')';
843 static auto serializationTestItems()
846 std::pair { Item::Oil , 876.54 },
847 std::pair { Item::Gas , 321.09 },
848 std::pair { Item::Water , 987.65 },
849 std::pair { Item::Mixture , 975.31 },
850 std::pair { Item::MixtureWithExponents, 765.43 },
862 SegmentPhaseQuantity velocity{};
863 SegmentPhaseQuantity holdup{};
864 SegmentPhaseQuantity viscosity{};
865 SegmentPhaseDensity density{};
866 std::size_t segNumber{};
868 bool operator==(
const Segment& seg2)
const
870 return (rates == seg2.rates)
871 && (pressures == seg2.pressures)
872 && (velocity == seg2.velocity)
873 && (holdup == seg2.holdup)
874 && (viscosity == seg2.viscosity)
875 && (density == seg2.density)
876 && (segNumber == seg2.segNumber);
879 template <
class MessageBufferType>
880 void write(MessageBufferType& buffer)
const;
882 template <
class MessageBufferType>
883 void read(MessageBufferType& buffer);
885 template <
class Serializer>
888 serializer(this->rates);
889 serializer(this->pressures);
890 serializer(this->velocity);
891 serializer(this->holdup);
892 serializer(this->viscosity);
893 serializer(this->density);
894 serializer(this->segNumber);
897 static Segment serializationTestObject()
900 Rates::serializationTestObject(),
901 SegmentPressures::serializationTestObject(),
902 SegmentPhaseQuantity::serializationTestObject(),
903 SegmentPhaseQuantity::serializationTestObject(),
904 SegmentPhaseQuantity::serializationTestObject(),
905 SegmentPhaseDensity::serializationTestObject(),
913 bool isProducer{
true};
915 ::Opm::WellProducerCMode prod {
916 ::Opm::WellProducerCMode::CMODE_UNDEFINED
919 ::Opm::WellInjectorCMode inj {
920 ::Opm::WellInjectorCMode::CMODE_UNDEFINED
925 return (this->isProducer == rhs.isProducer)
926 && ((this->isProducer && (this->prod == rhs.prod)) ||
927 (!this->isProducer && (this->inj == rhs.inj)));
930 template <
class MessageBufferType>
931 void write(MessageBufferType& buffer)
const;
933 template <
class MessageBufferType>
934 void read(MessageBufferType& buffer);
936 template<
class Serializer>
939 serializer(isProducer);
947 ::Opm::WellProducerCMode::BHP,
948 ::Opm::WellInjectorCMode::GRUP
956 enum class Quantity { WBP, WBP4, WBP5, WBP9 };
958 double& operator[](
const Quantity q)
960 return this->wbp_[
static_cast<std::size_t
>(q)];
963 double operator[](
const Quantity q)
const
965 return this->wbp_[
static_cast<std::size_t
>(q)];
970 return this->wbp_ == that.wbp_;
973 template <
class MessageBufferType>
974 void write(MessageBufferType& buffer)
const;
976 template <
class MessageBufferType>
977 void read(MessageBufferType& buffer);
979 template <
class Serializer>
982 serializer(this->wbp_);
989 wbp[Quantity::WBP] = 17.29;
990 wbp[Quantity::WBP4] = 2.718;
991 wbp[Quantity::WBP5] = 3.1415;
992 wbp[Quantity::WBP9] = 1.618;
998 static constexpr auto NumQuantities =
999 static_cast<std::size_t
>(Quantity::WBP9) + 1;
1001 std::array<double, NumQuantities> wbp_{};
1008 double concentration{0.};
1010 template<
class Serializer>
1014 serializer(concentration);
1018 return this->rate == filtrate.rate
1019 && this->total == filtrate.total
1020 && this->concentration == filtrate.concentration;
1027 res.concentration = 0.;
1031 template <
class MessageBufferType>
1032 void write(MessageBufferType& buffer)
const;
1034 template <
class MessageBufferType>
1035 void read(MessageBufferType& buffer);
1041 Bhp, OilRate, WaterRate, GasRate, ResVRate, LiquidRate,
1047 static std::string itemName(
const Item p)
1050 case Item::Bhp:
return "Bhp";
1051 case Item::OilRate:
return "OilRate";
1052 case Item::WaterRate:
return "WaterRate";
1053 case Item::GasRate:
return "GasRate";
1054 case Item::ResVRate:
return "ResVRate";
1055 case Item::LiquidRate:
return "LiquidRate";
1057 case Item::NumItems:
1058 return "Out of bounds (NumItems)";
1061 return "Unknown (" + std::to_string(
static_cast<int>(p)) +
')';
1064 static auto serializationTestItems()
1066 return std::vector {
1067 std::pair { Item::Bhp , 321.09 },
1068 std::pair { Item::OilRate , 987.65 },
1069 std::pair { Item::WaterRate , 975.31 },
1070 std::pair { Item::GasRate , 765.43 },
1071 std::pair { Item::ResVRate , 876.54 },
1072 std::pair { Item::LiquidRate, 54.32 },
1085 double temperature{0.0};
1087 double efficiency_scaling_factor{1.0};
1091 ::Opm::WellStatus dynamicStatus { Opm::WellStatus::OPEN };
1093 std::vector<Connection> connections{};
1094 std::unordered_map<std::size_t, Segment> segments{};
1097 WellControlLimits limits{};
1099 inline bool flowing()
const noexcept;
1101 template <
class MessageBufferType>
1102 void write(MessageBufferType& buffer)
const;
1104 template <
class MessageBufferType>
1105 void read(MessageBufferType& buffer);
1108 find_connection(
const Connection::global_index connection_grid_index)
const
1110 const auto connection =
1111 std::ranges::find_if(this->connections,
1113 {
return c.index == connection_grid_index; });
1115 if (connection == this->connections.end()) {
1119 return &*connection;
1123 find_connection(
const Connection::global_index connection_grid_index)
1125 const auto connection =
1126 std::ranges::find_if(this->connections,
1128 {
return c.index == connection_grid_index; });
1130 if (connection == this->connections.end()) {
1134 return &*connection;
1137 bool operator==(
const Well& well2)
const
1139 return (this->rates == well2.rates)
1140 && (this->bhp == well2.bhp)
1141 && (this->thp == well2.thp)
1142 && (this->temperature == well2.temperature)
1143 && (this->filtrate == well2.filtrate)
1144 && (this->control == well2.control)
1145 && (this->dynamicStatus == well2.dynamicStatus)
1146 && (this->connections == well2.connections)
1147 && (this->segments == well2.segments)
1148 && (this->current_control == well2.current_control)
1149 && (this->guide_rates == well2.guide_rates)
1150 && (this->limits == well2.limits)
1154 bool operator!=(
const Well& well2)
const
1156 return !(*
this == well2);
1159 template<
class Serializer>
1165 serializer(temperature);
1166 serializer(control);
1167 serializer(efficiency_scaling_factor);
1168 serializer(filtrate);
1169 serializer(dynamicStatus);
1170 serializer(connections);
1171 serializer(segments);
1172 serializer(current_control);
1173 serializer(guide_rates);
1177 static Well serializationTestObject()
1180 Rates::serializationTestObject(),
1186 WellFiltrate::serializationTestObject(),
1187 ::Opm::WellStatus::SHUT,
1188 {Connection::serializationTestObject()},
1189 {{0, Segment::serializationTestObject()}},
1190 CurrentControl::serializationTestObject(),
1191 GuideRateValue::serializationTestObject(),
1192 WellControlLimits::serializationTestObject()
1197 class Wells:
public std::map<std::string , Well> {
1200 double get(
const std::string& well_name , Rates::opt m)
const {
1201 const auto& well = this->find( well_name );
1202 if( well == this->end() )
return 0.0;
1204 return well->second.rates.get( m, 0.0 );
1207 double get(
const std::string& well_name , Rates::opt m,
const std::string& tracer_name)
const {
1208 const auto& well = this->find( well_name );
1209 if( well == this->end() )
return 0.0;
1211 return well->second.rates.get( m, 0.0, tracer_name);
1214 double get(
const std::string& well_name , Connection::global_index connection_grid_index, Rates::opt m)
const {
1215 const auto& witr = this->find( well_name );
1216 if( witr == this->end() )
return 0.0;
1218 const auto& well = witr->second;
1219 const auto connection =
1220 std::ranges::find_if(well.connections,
1222 { return c.index == connection_grid_index; });
1224 if (connection == well.connections.end()) {
1228 return connection->rates.get(m, 0.0);
1231 template <
class MessageBufferType>
1232 void write(MessageBufferType& buffer)
const {
1233 unsigned int size = this->size();
1235 for (
const auto& witr : *
this) {
1236 const std::string& name = witr.first;
1238 const Well& well = witr.second;
1243 template <
class MessageBufferType>
1244 void read(MessageBufferType& buffer) {
1247 for (
size_t i = 0; i < size; ++i) {
1252 auto result = this->emplace(name, well);
1255 if (!result.second && result.first->second != well) {
1256 OPM_THROW(std::runtime_error,
"Received different output data for well " + name +
" from more than one process, the output of this simulation will be wrong!");
1257 }
else if (!result.second) {
1258 OpmLog::warning(
"Received consistently duplicated output data for well " + name +
" from more than one process - this might be problematic!");
1263 template<
class Serializer>
1266 serializer(
static_cast<std::map<std::string,Well>&
>(*
this));
1269 static Wells serializationTestObject()
1272 w.insert({
"test_well", Well::serializationTestObject()});
1280 std::unordered_map<std::string, WellBlockAvgPress> values{};
1282 template <
class MessageBufferType>
1283 void write(MessageBufferType& buffer)
const;
1285 template <
class MessageBufferType>
1286 void read(MessageBufferType& buffer);
1290 return this->values == that.values;
1293 template <
class Serializer>
1296 serializer(this->values);
1302 { {
"I-45", WellBlockAvgPress::serializationTestObject() } },
1310 const auto mand =
static_cast< enum_size
>( this->mask )
1311 &
static_cast< enum_size
>( m );
1313 return static_cast< opt
>( mand ) == m;
1317 if( !this->
has( m ) )
1318 throw std::invalid_argument(
"Uninitialized value." );
1320 return this->get_ref( m );
1324 if( !this->
has( m ) )
return default_value;
1326 return this->get_ref( m );
1329 inline double Rates::get( opt m,
double default_value,
const std::string& tracer_name)
const {
1330 if( !this->
has( m ) )
return default_value;
1332 if( m == opt::tracer && this->tracer.find(tracer_name) == this->tracer.end())
return default_value;
1334 return this->get_ref( m, tracer_name);
1338 this->get_ref( m ) = value;
1341 this->mask =
static_cast< opt
>(
1342 static_cast< enum_size
>( this->mask ) |
1343 static_cast< enum_size
>( m )
1349 inline Rates&
Rates::set( opt m,
double value ,
const std::string& tracer_name ) {
1350 this->get_ref( m , tracer_name) = value;
1353 this->mask =
static_cast< opt
>(
1354 static_cast< enum_size
>( this->mask ) |
1355 static_cast< enum_size
>( m )
1361 inline bool Rates::operator==(
const Rates& rate)
const
1363 return mask == rate.mask &&
1367 polymer == rate.polymer &&
1368 solvent == rate.solvent &&
1369 energy == rate.energy &&
1370 dissolved_gas == rate.dissolved_gas &&
1371 vaporized_oil == rate.vaporized_oil &&
1372 reservoir_water == rate.reservoir_water &&
1373 reservoir_oil == rate.reservoir_oil &&
1374 reservoir_gas == rate.reservoir_gas &&
1375 productivity_index_water == rate.productivity_index_water &&
1376 productivity_index_gas == rate.productivity_index_gas &&
1377 productivity_index_oil == rate.productivity_index_oil &&
1378 well_potential_water == rate.well_potential_water &&
1379 well_potential_oil == rate.well_potential_oil &&
1380 well_potential_gas == rate.well_potential_gas &&
1381 brine == rate.brine &&
1383 tracer == rate.tracer &&
1384 microbial == rate.microbial &&
1385 oxygen == rate.oxygen &&
1386 urea == rate.urea &&
1387 vaporized_water == rate.vaporized_water &&
1388 mass_gas == rate.mass_gas &&
1389 mass_wat == rate.mass_wat &&
1390 wat_frac == rate.wat_frac;
1403 inline const double& Rates::get_ref( opt m )
const {
1405 case opt::wat:
return this->wat;
1406 case opt::oil:
return this->oil;
1407 case opt::gas:
return this->gas;
1408 case opt::polymer:
return this->polymer;
1409 case opt::solvent:
return this->solvent;
1410 case opt::energy:
return this->energy;
1411 case opt::dissolved_gas:
return this->dissolved_gas;
1412 case opt::vaporized_oil:
return this->vaporized_oil;
1413 case opt::reservoir_water:
return this->reservoir_water;
1414 case opt::reservoir_oil:
return this->reservoir_oil;
1415 case opt::reservoir_gas:
return this->reservoir_gas;
1416 case opt::productivity_index_water:
return this->productivity_index_water;
1417 case opt::productivity_index_oil:
return this->productivity_index_oil;
1418 case opt::productivity_index_gas:
return this->productivity_index_gas;
1419 case opt::well_potential_water:
return this->well_potential_water;
1420 case opt::well_potential_oil:
return this->well_potential_oil;
1421 case opt::well_potential_gas:
return this->well_potential_gas;
1422 case opt::brine:
return this->brine;
1423 case opt::alq:
return this->alq;
1426 case opt::microbial:
return this->microbial;
1427 case opt::oxygen:
return this->oxygen;
1428 case opt::urea:
return this->urea;
1429 case opt::vaporized_water:
return this->vaporized_water;
1430 case opt::mass_gas:
return this->mass_gas;
1431 case opt::mass_wat:
return this->mass_wat;
1432 case opt::wat_frac:
return this->wat_frac;
1435 throw std::invalid_argument(
1436 "Unknown value type '"
1437 + std::to_string(
static_cast< enum_size
>( m ) )
1442 inline const double& Rates::get_ref( opt m,
const std::string& tracer_name )
const {
1443 if (m != opt::tracer)
1444 throw std::logic_error(
"Logic error - should be called with tracer argument");
1446 return this->tracer.at(tracer_name);
1449 inline double& Rates::get_ref( opt m ) {
1450 return const_cast< double&
>(
1451 static_cast< const Rates*
>( this )->get_ref( m )
1455 inline double& Rates::get_ref( opt m,
const std::string& tracer_name ) {
1456 if (m == opt::tracer) this->tracer.emplace(tracer_name, 0.0);
1457 return this->tracer.at(tracer_name);
1461 return ((this->wat != 0) ||
1466 inline bool Well::flowing() const noexcept {
1470 template <
class MessageBufferType>
1471 void Rates::write(MessageBufferType& buffer)
const {
1472 buffer.write(this->mask);
1473 buffer.write(this->wat);
1474 buffer.write(this->oil);
1475 buffer.write(this->gas);
1476 buffer.write(this->polymer);
1477 buffer.write(this->solvent);
1478 buffer.write(this->energy);
1479 buffer.write(this->dissolved_gas);
1480 buffer.write(this->vaporized_oil);
1481 buffer.write(this->reservoir_water);
1482 buffer.write(this->reservoir_oil);
1483 buffer.write(this->reservoir_gas);
1484 buffer.write(this->productivity_index_water);
1485 buffer.write(this->productivity_index_oil);
1486 buffer.write(this->productivity_index_gas);
1487 buffer.write(this->well_potential_water);
1488 buffer.write(this->well_potential_oil);
1489 buffer.write(this->well_potential_gas);
1490 buffer.write(this->brine);
1491 buffer.write(this->alq);
1494 unsigned int size = this->tracer.size();
1496 for (
const auto& [name, rate] : this->tracer) {
1501 buffer.write(this->microbial);
1502 buffer.write(this->oxygen);
1503 buffer.write(this->urea);
1504 buffer.write(this->vaporized_water);
1505 buffer.write(this->mass_gas);
1506 buffer.write(this->mass_wat);
1507 buffer.write(this->wat_frac);
1510 template <
class MessageBufferType>
1511 void ConnectionFiltrate::write(MessageBufferType& buffer)
const {
1512 buffer.write(this->rate);
1513 buffer.write(this->total);
1514 buffer.write(this->skin_factor);
1515 buffer.write(this->thickness);
1516 buffer.write(this->perm);
1517 buffer.write(this->poro);
1518 buffer.write(this->radius);
1519 buffer.write(this->area_of_flow);
1520 buffer.write(this->flow_factor);
1521 buffer.write(this->fracture_rate);
1524 template <
class MessageBufferType>
1525 void Connection::write(MessageBufferType& buffer)
const {
1526 buffer.write(this->index);
1527 this->rates.write(buffer);
1528 buffer.write(this->pressure);
1529 buffer.write(this->reservoir_rate);
1530 buffer.write(this->cell_pressure);
1531 buffer.write(this->cell_saturation_water);
1532 buffer.write(this->cell_saturation_gas);
1533 buffer.write(this->effective_Kh);
1534 buffer.write(this->trans_factor);
1535 buffer.write(this->d_factor);
1536 buffer.write(this->compact_mult);
1537 buffer.write(this->lgr_grid);
1538 this->filtrate.write(buffer);
1539 this->fracture.write(buffer);
1540 this->
fract.write(buffer);
1543 template <
class MessageBufferType>
1544 void Segment::write(MessageBufferType& buffer)
const
1546 buffer.write(this->segNumber);
1547 this->rates.write(buffer);
1548 this->pressures.write(buffer);
1549 this->velocity.write(buffer);
1550 this->holdup.write(buffer);
1551 this->viscosity.write(buffer);
1552 this->density.write(buffer);
1555 template <
class MessageBufferType>
1556 void CurrentControl::write(MessageBufferType& buffer)
const
1558 buffer.write(this->isProducer);
1559 if (this->isProducer) {
1560 buffer.write(this->prod);
1563 buffer.write(this->inj);
1567 template <
class MessageBufferType>
1568 void WellBlockAvgPress::write(MessageBufferType& buffer)
const
1570 for (
const auto& quantity : this->wbp_) {
1571 buffer.write(quantity);
1575 template <
class MessageBufferType>
1576 void WellFiltrate::write(MessageBufferType& buffer)
const
1578 buffer.write(this->rate);
1579 buffer.write(this->total);
1580 buffer.write(this->concentration);
1583 template <
class MessageBufferType>
1584 void Well::write(MessageBufferType& buffer)
const
1586 this->rates.write(buffer);
1588 buffer.write(this->bhp);
1589 buffer.write(this->thp);
1590 buffer.write(this->temperature);
1591 buffer.write(this->control);
1592 buffer.write(this->efficiency_scaling_factor);
1594 this->filtrate.write(buffer);
1597 const auto status = ::Opm::WellStatus2String(this->dynamicStatus);
1598 buffer.write(status);
1602 const unsigned int size = this->connections.size();
1605 for (
const Connection& comp : this->connections) {
1612 static_cast<unsigned int>(this->segments.size());
1615 for (
const auto& seg : this->segments) {
1616 seg.second.write(buffer);
1620 this->current_control.write(buffer);
1621 this->guide_rates.write(buffer);
1622 this->limits.write(buffer);
1625 template <
class MessageBufferType>
1626 void WellBlockAveragePressures::write(MessageBufferType& buffer)
const
1628 buffer.write(this->values.size());
1630 for (
const auto& [well, value] : this->values) {
1632 value.write(buffer);
1636 template <
class MessageBufferType>
1637 void Rates::read(MessageBufferType& buffer) {
1638 buffer.read(this->mask);
1639 buffer.read(this->wat);
1640 buffer.read(this->oil);
1641 buffer.read(this->gas);
1642 buffer.read(this->polymer);
1643 buffer.read(this->solvent);
1644 buffer.read(this->energy);
1645 buffer.read(this->dissolved_gas);
1646 buffer.read(this->vaporized_oil);
1647 buffer.read(this->reservoir_water);
1648 buffer.read(this->reservoir_oil);
1649 buffer.read(this->reservoir_gas);
1650 buffer.read(this->productivity_index_water);
1651 buffer.read(this->productivity_index_oil);
1652 buffer.read(this->productivity_index_gas);
1653 buffer.read(this->well_potential_water);
1654 buffer.read(this->well_potential_oil);
1655 buffer.read(this->well_potential_gas);
1656 buffer.read(this->brine);
1657 buffer.read(this->alq);
1662 for (
size_t i = 0; i < size; ++i) {
1663 std::string tracer_name;
1664 buffer.read(tracer_name);
1666 buffer.read(tracer_rate);
1667 this->tracer.emplace(tracer_name, tracer_rate);
1670 buffer.read(this->microbial);
1671 buffer.read(this->oxygen);
1672 buffer.read(this->urea);
1673 buffer.read(this->vaporized_water);
1674 buffer.read(this->mass_gas);
1675 buffer.read(this->mass_wat);
1676 buffer.read(this->wat_frac);
1679 template <
class MessageBufferType>
1680 void ConnectionFiltrate::read(MessageBufferType& buffer) {
1681 buffer.read(this->rate);
1682 buffer.read(this->total);
1683 buffer.read(this->skin_factor);
1684 buffer.read(this->thickness);
1685 buffer.read(this->perm);
1686 buffer.read(this->poro);
1687 buffer.read(this->radius);
1688 buffer.read(this->area_of_flow);
1689 buffer.read(this->flow_factor);
1690 buffer.read(this->fracture_rate);
1693 template <
class MessageBufferType>
1694 void Connection::read(MessageBufferType& buffer) {
1695 buffer.read(this->index);
1696 this->rates.read(buffer);
1697 buffer.read(this->pressure);
1698 buffer.read(this->reservoir_rate);
1699 buffer.read(this->cell_pressure);
1700 buffer.read(this->cell_saturation_water);
1701 buffer.read(this->cell_saturation_gas);
1702 buffer.read(this->effective_Kh);
1703 buffer.read(this->trans_factor);
1704 buffer.read(this->d_factor);
1705 buffer.read(this->compact_mult);
1706 buffer.read(this->lgr_grid);
1707 this->filtrate.read(buffer);
1708 this->fracture.read(buffer);
1709 this->
fract.read(buffer);
1712 template <
class MessageBufferType>
1713 void Segment::read(MessageBufferType& buffer)
1715 buffer.read(this->segNumber);
1716 this->rates.read(buffer);
1717 this->pressures.read(buffer);
1718 this->velocity.read(buffer);
1719 this->holdup.read(buffer);
1720 this->viscosity.read(buffer);
1721 this->density.read(buffer);
1724 template <
class MessageBufferType>
1725 void CurrentControl::read(MessageBufferType& buffer)
1727 buffer.read(this->isProducer);
1728 if (this->isProducer) {
1729 buffer.read(this->prod);
1732 buffer.read(this->inj);
1736 template <
class MessageBufferType>
1737 void WellBlockAvgPress::read(MessageBufferType& buffer)
1739 for (
auto& quantity : this->wbp_) {
1740 buffer.read(quantity);
1744 template <
class MessageBufferType>
1745 void WellFiltrate::read(MessageBufferType& buffer)
1747 buffer.read(this->rate);
1748 buffer.read(this->total);
1749 buffer.read(this->concentration);
1752 template <
class MessageBufferType>
1753 void Well::read(MessageBufferType& buffer)
1755 this->rates.read(buffer);
1757 buffer.read(this->bhp);
1758 buffer.read(this->thp);
1759 buffer.read(this->temperature);
1760 buffer.read(this->control);
1761 buffer.read(this->efficiency_scaling_factor);
1763 this->filtrate.read(buffer);
1766 auto status = std::string{};
1767 buffer.read(status);
1768 this->dynamicStatus = ::Opm::WellStatusFromString(status);
1773 unsigned int size = 0;
1776 this->connections.resize(size);
1777 for (
auto& connection : this->connections) {
1778 connection.read(buffer);
1783 const auto nSeg = [&buffer]() ->
unsigned int
1791 for (
auto segID = 0*nSeg; segID < nSeg; ++segID) {
1792 auto seg = Segment{};
1795 const auto segNumber = seg.segNumber;
1796 this->segments.emplace(segNumber, std::move(seg));
1799 this->current_control.read(buffer);
1800 this->guide_rates.read(buffer);
1801 this->limits.read(buffer);
1804 template <
class MessageBufferType>
1805 void WellBlockAveragePressures::read(MessageBufferType& buffer)
1807 const auto numWells = [&buffer,
this]()
1809 auto size = 0*this->values.size();
1815 auto wellName = std::string{};
1816 for (
auto well = 0*numWells; well < numWells; ++well) {
1817 buffer.read(wellName);
1819 this->values[wellName].read(buffer);
Class for (de-)serializing.
Definition Serializer.hpp:94
Definition GuideRateValue.hpp:31
bool flowing() const
Returns true if any of the rates oil, gas, water is nonzero.
Definition Wells.hpp:1460
double get(opt m) const
Read the value indicated by m.
Definition Wells.hpp:1316
Rates & set(opt m, double value)
Set the value specified by m.
Definition Wells.hpp:1337
bool has(opt) const
Query if a value is set.
Definition Wells.hpp:1309
Definition Wells.hpp:1197
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Statistics collection for a single quantity.
Definition Wells.hpp:365
void serializeOp(Serializer &serializer)
Convert between byte array and object representation.
Definition Wells.hpp:394
void write(MessageBufferType &buffer) const
MPI communication protocol–serialisation operation.
Definition Wells.hpp:420
double stdev
Unbiased sample standard deviation.
Definition Wells.hpp:378
double avg
Arithmetic average.
Definition Wells.hpp:367
double min
Minimum value.
Definition Wells.hpp:373
void read(MessageBufferType &buffer)
MPI communication protocol–deserialisation operation.
Definition Wells.hpp:430
double max
Maximum value.
Definition Wells.hpp:370
bool operator==(const Statistics &that) const
Equality predicate.
Definition Wells.hpp:409
static Statistics serializationTestObject()
Create a serialization test object.
Definition Wells.hpp:381
Connection Level Fracturing Statistics.
Definition Wells.hpp:362
void read(MessageBufferType &buffer)
MPI communication protocol–deserialisation operation.
Definition Wells.hpp:508
void write(MessageBufferType &buffer) const
MPI communication protocol–serialisation operation.
Definition Wells.hpp:498
static ConnectionFracturing serializationTestObject()
Create a serialisation test object.
Definition Wells.hpp:454
std::size_t numCells
Sample size.
Definition Wells.hpp:442
bool operator==(const ConnectionFracturing &that) const
Equality predicate.
Definition Wells.hpp:487
Statistics width
Statistical measures for connection's fracture fracture width.
Definition Wells.hpp:451
void serializeOp(Serializer &serializer)
Convert between byte array and object representation.
Definition Wells.hpp:472
Statistics rate
Statistical measures for connection's fracture fracture flow rate.
Definition Wells.hpp:448
Statistics press
Statistical measures for connection's fracture pressures.
Definition Wells.hpp:445
ConnectionFracturing fract
Connection level fracturing statistics.
Definition Wells.hpp:541
Definition Wells.hpp:1279
Definition Wells.hpp:1039
Definition Wells.hpp:1005
Definition Wells.hpp:1080