/* Copyright 2014 SINTEF ICT, Applied Mathematics. Copyright 2017 IRIS AS This file is part of the Open Porous Media project (OPM). OPM 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. OPM 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 OPM. If not, see . */ #ifndef OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED #define OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED #include #include #include #include #include #include #include #include #include #include #include #include #include namespace Opm { class ParallelWellInfo; class Schedule; /// The state of a set of wells, tailored for use by the fully /// implicit blackoil simulator. class WellStateFullyImplicitBlackoil : public WellState { typedef WellState BaseType; public: static const uint64_t event_mask = ScheduleEvents::WELL_STATUS_CHANGE + ScheduleEvents::PRODUCTION_UPDATE + ScheduleEvents::INJECTION_UPDATE; typedef BaseType :: WellMapType WellMapType; virtual ~WellStateFullyImplicitBlackoil() = default; // TODO: same definition with WellInterface, eventually they should go to a common header file. static const int Water = BlackoilPhases::Aqua; static const int Oil = BlackoilPhases::Liquid; static const int Gas = BlackoilPhases::Vapour; using BaseType :: wellRates; using BaseType :: bhp; using BaseType :: perfPress; using BaseType :: numPhases; using BaseType :: updateStatus; explicit WellStateFullyImplicitBlackoil(const PhaseUsage& pu) : WellState(pu) { } const WellMapType& wellMap() const { return wellMap_; } WellMapType& wellMap() { return wellMap_; } int numWells() const { return wellMap_.size(); } const ParallelWellInfo& parallelWellInfo(std::size_t well_index) const; /// Allocate and initialize if wells is non-null. Also tries /// to give useful initial values to the bhp(), wellRates() /// and perfPhaseRates() fields, depending on controls void init(const std::vector& cellPressures, const Schedule& schedule, const std::vector& wells_ecl, const std::vector& parallel_well_info, const int report_step, const WellStateFullyImplicitBlackoil* prevState, const std::vector>& well_perf_data, const SummaryState& summary_state); void resize(const std::vector& wells_ecl, const std::vector& parallel_well_info, const Schedule& schedule, const bool handle_ms_well, const size_t numCells, const std::vector>& well_perf_data, const SummaryState& summary_state); /// One rate per phase and well connection. std::vector& mutable_perfPhaseRates() { return perfphaserates_; } const std::vector& perfPhaseRates() const { return perfphaserates_; } /// One current control per injecting well. Well::InjectorCMode currentInjectionControl(std::size_t well_index) const { return current_injection_controls_[well_index]; } void currentInjectionControl(std::size_t well_index, Well::InjectorCMode cmode) { current_injection_controls_[well_index] = cmode; } /// One current control per producing well. Well::ProducerCMode currentProductionControl(std::size_t well_index) const { return current_production_controls_[well_index]; } void currentProductionControl(std::size_t well_index, Well::ProducerCMode cmode) { current_production_controls_[well_index] = cmode; } void setCurrentWellRates(const std::string& wellName, const std::vector& rates ) { well_rates[wellName].second = rates; } const std::vector& currentWellRates(const std::string& wellName) const; bool hasWellRates(const std::string& wellName) const { return this->well_rates.find(wellName) != this->well_rates.end(); } template void gatherVectorsOnRoot(const std::vector< data::Connection >& from_connections, std::vector< data::Connection >& to_connections, const Communication& comm) const; data::Wells report(const int* globalCellIdxMap, const std::function& wasDynamicallyClosed) const; void reportConnections(data::Well& well, const PhaseUsage &pu, const WellMapType::value_type& wt, const int* globalCellIdxMap) const; /// init the MS well related. void initWellStateMSWell(const std::vector& wells_ecl, const WellStateFullyImplicitBlackoil* prev_well_state); static void calculateSegmentRates(const std::vector>& segment_inlets, const std::vector>&segment_perforations, const std::vector& perforation_rates, const int np, const int segment, std::vector& segment_rates); Events& events(std::size_t well_index) { return this->events_[well_index]; } const std::vector& firstPerfIndex() const { return first_perf_index_; } /// One rate pr well connection. std::vector& perfRateSolvent() { return perfRateSolvent_; } const std::vector& perfRateSolvent() const { return perfRateSolvent_; } /// One rate pr well double solventWellRate(const int w) const; /// One rate pr well connection. std::vector& perfRatePolymer() { return perfRatePolymer_; } const std::vector& perfRatePolymer() const { return perfRatePolymer_; } /// One rate pr well double polymerWellRate(const int w) const; /// One rate pr well connection. std::vector& perfRateBrine() { return perfRateBrine_; } const std::vector& perfRateBrine() const { return perfRateBrine_; } /// One rate pr well double brineWellRate(const int w) const; const WellContainer>& wellReservoirRates() const { return well_reservoir_rates_; } std::vector& wellReservoirRates(std::size_t well_index) { return well_reservoir_rates_[well_index]; } const std::vector& wellReservoirRates(std::size_t well_index) const { return well_reservoir_rates_[well_index]; } std::vector& wellDissolvedGasRates() { return well_dissolved_gas_rates_; } std::vector& wellVaporizedOilRates() { return well_vaporized_oil_rates_; } const std::vector& segRates() const { return seg_rates_; } std::vector& segRates() { return seg_rates_; } const std::vector& segPress() const { return seg_press_; } std::vector& segPressDrop() { return seg_pressdrop_; } const std::vector& segPressDrop() const { return seg_pressdrop_; } std::vector& segPressDropFriction() { return seg_pressdrop_friction_; } const std::vector& segPressDropFriction() const { return seg_pressdrop_friction_; } std::vector& segPressDropHydroStatic() { return seg_pressdrop_hydorstatic_; } const std::vector& segPressDropHydroStatic() const { return seg_pressdrop_hydorstatic_; } std::vector& segPressDropAcceleration() { return seg_pressdrop_acceleration_; } const std::vector& segPressDropAcceleration() const { return seg_pressdrop_acceleration_; } std::vector& segPress() { return seg_press_; } int numSegment() const { return nseg_; } int topSegmentIndex(const int w) const; std::vector& productivityIndex() { return productivity_index_; } const std::vector& productivityIndex() const { return productivity_index_; } std::vector& connectionProductivityIndex() { return this->conn_productivity_index_; } const std::vector& connectionProductivityIndex() const { return this->conn_productivity_index_; } std::vector& wellPotentials() { return well_potentials_; } const std::vector& wellPotentials() const { return well_potentials_; } std::vector& perfThroughput() { return perf_water_throughput_; } const std::vector& perfThroughput() const { return perf_water_throughput_; } std::vector& perfSkinPressure() { return perf_skin_pressure_; } const std::vector& perfSkinPressure() const { return perf_skin_pressure_; } std::vector& perfWaterVelocity() { return perf_water_velocity_; } const std::vector& perfWaterVelocity() const { return perf_water_velocity_; } void shutWell(int well_index) override; template void communicateGroupRates(const Comm& comm); template void updateGlobalIsGrup(const Comm& comm); bool isInjectionGrup(const std::string& name) const { return this->global_well_info.value().in_injecting_group(name); } bool isProductionGrup(const std::string& name) const { return this->global_well_info.value().in_producing_group(name); } double getALQ( const std::string& name) const { return this->alq_state.get(name); } void setALQ( const std::string& name, double value) { this->alq_state.set(name, value); } bool gliftCheckAlqOscillation(const std::string &name) const { return this->alq_state.oscillation(name); } int gliftGetAlqDecreaseCount(const std::string &name) { return this->alq_state.get_decrement_count(name); } int gliftGetAlqIncreaseCount(const std::string &name) { return this->alq_state.get_increment_count(name); } void gliftUpdateAlqIncreaseCount(const std::string &name, bool increase) { this->alq_state.update_count(name, increase); } bool gliftOptimizationEnabled() const { return do_glift_optimization_; } void gliftTimeStepInit() { this->alq_state.reset_count(); disableGliftOptimization(); } void disableGliftOptimization() { do_glift_optimization_ = false; } void enableGliftOptimization() { do_glift_optimization_ = true; } int wellNameToGlobalIdx(const std::string &name) { return this->global_well_info.value().well_index(name); } std::string globalIdxToWellName(const int index) { return this->global_well_info.value().well_name(index); } bool wellIsOwned(std::size_t well_index, const std::string& wellName) const; bool wellIsOwned(const std::string& wellName) const; /// Special purpose method to support dynamically rescaling a well's /// CTFs through WELPI. /// /// \param[in] well_index Process-local linear index of single well. /// Must be in the range 0..numWells()-1. /// /// \param[in] well_perf_data New perforation data. Only /// PerforationData::connection_transmissibility_factor actually /// used (overwrites existing internal values). void resetConnectionTransFactors(const int well_index, const std::vector& well_perf_data); private: std::vector perfphaserates_; WellContainer is_producer_; // Size equal to number of local wells. // vector with size number of wells +1. // iterate over all perforations of a given well // for (int perf = first_perf_index_[well_index]; perf < first_perf_index_[well_index] + num_perf_[well_index]; ++perf) std::vector first_perf_index_; std::vector num_perf_; WellContainer current_injection_controls_; WellContainer current_production_controls_; // Use of std::optional<> here is a technical crutch, the // WellStateFullyImplicitBlackoil class should be default constructible, // whereas the GlobalWellInfo is not. std::optional global_well_info; std::map>> well_rates; ALQState alq_state; bool do_glift_optimization_; std::vector perfRateSolvent_; // only for output std::vector perfRatePolymer_; std::vector perfRateBrine_; // it is the throughput of water flow through the perforations // it is used as a measure of formation damage around well-bore due to particle deposition // it will only be used for injectors to check the injectivity std::vector perf_water_throughput_; // skin pressure of peforation // it will only be used for injectors to check the injectivity std::vector perf_skin_pressure_; // it will only be used for injectors to check the injectivity // water velocity of perforation std::vector perf_water_velocity_; // phase rates under reservoir condition for wells // or voidage phase rates WellContainer> well_reservoir_rates_; // dissolved gas rates or solution gas production rates // should be zero for injection wells std::vector well_dissolved_gas_rates_; // vaporized oil rates or solution oil producation rates // should be zero for injection wells std::vector well_vaporized_oil_rates_; // some events happens to the well, like this well is a new well // or new well control keywords happens // \Note: for now, only WCON* keywords, and well status change is considered WellContainer events_; // MS well related // for StandardWell, the number of segments will be one std::vector seg_rates_; std::vector seg_press_; // the index of the top segments, which is used to locate the // multisegment well related information in WellState std::vector top_segment_index_; int nseg_; // total number of the segments // The following data are only recorded for output // pressure drop std::vector seg_pressdrop_; // frictional pressure drop std::vector seg_pressdrop_friction_; // hydrostatic pressure drop std::vector seg_pressdrop_hydorstatic_; // accelerational pressure drop std::vector seg_pressdrop_acceleration_; // Productivity Index std::vector productivity_index_; // Connection-level Productivity Index std::vector conn_productivity_index_; // Well potentials std::vector well_potentials_; /// Map segment index to segment number, mostly for MS wells. /// /// Segment number (one-based) of j-th segment in i-th well is /// \code /// const auto top = topSegmentIndex(i); /// const auto seg_No = seg_number_[top + j]; /// \end std::vector seg_number_; data::Segment reportSegmentResults(const PhaseUsage& pu, const int well_id, const int seg_ix, const int seg_no) const; int numSegments(const int well_id) const; int segmentNumber(const int well_id, const int seg_id) const; // If the ALQ has changed since the previous report step, // reset current_alq and update default_alq. ALQ is used for // constant lift gas injection and for gas lift optimization // (THP controlled wells). // // NOTE: If a well is no longer used (e.g. it is shut down) // it is still kept in the maps "default_alq_" and "current_alq_". Since the // number of unused entries should be small (negligible memory // overhead) this is simpler than writing code to delete it. // void updateWellsDefaultALQ(const std::vector& wells_ecl); }; } // namespace Opm #endif // OPM_WELLSTATEFULLYIMPLICITBLACKOIL_HEADER_INCLUDED