/* Copyright 2016 SINTEF ICT, Applied Mathematics. Copyright 2016 Statoil ASA. 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_MULTISEGMENTWELLS_HEADER_INCLUDED #define OPM_MULTISEGMENTWELLS_HEADER_INCLUDED #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace Opm { /// Class for handling the multi-segment well model class MultisegmentWells { public: // --------- Types --------- using ADB = AutoDiffBlock; using Vector = ADB::V; // Well operations and data needed. struct MultisegmentWellOps { explicit MultisegmentWellOps(const std::vector& wells_ms); Eigen::SparseMatrix w2p; // well -> perf (scatter) Eigen::SparseMatrix p2w; // perf -> well (gather) Eigen::SparseMatrix w2s; // well -> segment (scatter) Eigen::SparseMatrix s2w; // segment -> well (gather) Eigen::SparseMatrix s2p; // segment -> perf (scatter) Eigen::SparseMatrix p2s; // perf -> segment (gather) Eigen::SparseMatrix s2s_inlets; // segment -> its inlet segments Eigen::SparseMatrix s2s_outlet; // segment -> its outlet segment Eigen::SparseMatrix topseg2w; // top segment -> well AutoDiffMatrix eliminate_topseg; // change the top segment related to be zero std::vector well_cells; // the set of perforated cells Vector conn_trans_factors; // connection transmissibility factors bool has_multisegment_wells; // flag indicating whether there is any muli-segment well }; // copied from BlackoilModelBase // should put to somewhere better using DataBlock = Eigen::Array; // --------- Public methods --------- // TODO: using a vector of WellMultiSegmentConstPtr for now // TODO: it should use const Wells or something else later. MultisegmentWells(const Wells* wells_arg, const std::vector& wells_ecl, const int time_step); std::vector createMSWellVector(const Wells* wells_arg, const std::vector& wells_ecl, const int time_step); void init(const BlackoilPropsAdInterface* fluid_arg, const std::vector* active_arg, const std::vector* pc_arg, const VFPProperties* vfp_properties_arg, const double gravity_arg, const Vector& depth_arg); const std::vector& wells() const; const MultisegmentWellOps& wellOps() const; const Wells& wellsStruct() const; int numPhases() const { return num_phases_; }; int numWells() const { return wells_multisegment_.size(); } int numSegment() const { return nseg_total_; }; int numPerf() const { return nperf_total_; }; bool wellsActive() const { return wells_active_; }; void setWellsActive(const bool wells_active) { wells_active_ = wells_active; }; bool localWellsActive() const { return ! wells_multisegment_.empty(); } template void extractWellPerfProperties(const SolutionState& state, const std::vector& rq, std::vector& mob_perfcells, std::vector& b_perfcells) const; Vector& wellPerforationCellPressureDiffs() { return well_perforation_cell_pressure_diffs_; }; Vector& wellSegmentPerforationDepthDiffs() { return well_segment_perforation_depth_diffs_; }; const Vector& wellPerforationCellDensities() const { return well_perforation_cell_densities_; }; Vector& wellPerforationCellDensities() { return well_perforation_cell_densities_; }; const std::vector& segmentCompSurfVolumeInitial() const { return segment_comp_surf_volume_initial_; }; std::vector& segmentCompSurfVolumeInitial() { return segment_comp_surf_volume_initial_; }; const std::vector& segmentCompSurfVolumeCurrent() const { return segment_comp_surf_volume_current_; }; const std::vector& topWellSegments() const { return top_well_segments_; }; std::vector& topWellSegments() { return top_well_segments_; }; Vector& segVDt() { return segvdt_; }; const Vector& wellPerforationDensities() const { return well_perforation_densities_; }; Vector& wellPerforationDensities() { return well_perforation_densities_; }; const Vector& wellPerforationPressureDiffs() const { return well_perforation_pressure_diffs_; }; Vector& wellPerforationPressureDiffs() { return well_perforation_pressure_diffs_; }; template void updateWellState(const Vector& dwells, const double dpmaxrel, WellState& well_state) const; // TODO: some arguments can be removed later // TODO: compi will be required in the multisegment wells template void computeWellFlux(const SolutionState& state, const std::vector& mob_perfcells, const std::vector& b_perfcells, Vector& aliveWells, std::vector& cq_s) const; template void updatePerfPhaseRatesAndPressures(const std::vector& cq_s, const SolutionState& state, WellState& xw) const; // Calculate the density of the mixture in the segments // And the surface volume of the components in the segments by dt template void computeSegmentFluidProperties(const SolutionState& state); void computeSegmentPressuresDelta(const double grav); template void addWellFluxEq(const std::vector& cq_s, const SolutionState& state, LinearisedBlackoilResidual& residual); template void addWellControlEq(const SolutionState& state, const WellState& xw, const Vector& aliveWells, LinearisedBlackoilResidual& residual); template void updateWellControls(const bool terminal_output, WellState& xw) const; // TODO: these code are same with the StandardWells // to find a better solution later. void variableStateWellIndices(std::vector& indices, int& next) const; std::vector variableWellStateIndices() const; template void variableWellStateInitials(const WellState& xw, std::vector& vars0) const; template void computeWellConnectionPressures(const SolutionState& state, const WellState& xw, const std::vector& kr_adb, const std::vector& fluid_density); protected: // TODO: probably a wells_active_ will be required here. bool wells_active_; std::vector wells_multisegment_; MultisegmentWellOps wops_ms_; const int num_phases_; int nseg_total_; int nperf_total_; // TODO: put the Wells here to simplify the interface // TODO: at the moment, both wells_ and wells_multisegment_ // TODO: acutally contain all the wells // TODO: they should be split eventually. const Wells* wells_; const BlackoilPropsAdInterface* fluid_; const std::vector* active_; const std::vector* phase_condition_; const VFPProperties* vfp_properties_; double gravity_; // The depth of the all the cell centers // It is different from the perforation depth in MultisegmentWells Vector perf_cell_depth_; // Pressure correction due to the different depth of the perforation // and the cell center of the grid block // For the non-segmented wells, since the perforation are forced to be // at the center of the grid cell, it should be ZERO. // It only applies to the mutli-segmented wells. Vector well_perforation_cell_pressure_diffs_; // Pressure correction due to the depth differennce between segment depth and perforation depth. ADB well_segment_perforation_pressure_diffs_; // The depth difference between segment nodes and perforations Vector well_segment_perforation_depth_diffs_; // the average of the fluid densities in the grid block // which is used to calculate the hydrostatic head correction due to the depth difference of the perforation // and the cell center of the grid block Vector well_perforation_cell_densities_; // the density of the fluid mixture in the segments // which is calculated in an implicit way ADB well_segment_densities_; // the hydrostatic pressure drop between segment nodes // calculated with the above density of fluid mixtures // for the top segment, they should always be zero for the moment. ADB well_segment_pressures_delta_; // the surface volume of components in the segments // the initial value at the beginning of the time step std::vector segment_comp_surf_volume_initial_; // the value within the current iteration. std::vector segment_comp_surf_volume_current_; // the mass flow rate in the segments ADB segment_mass_flow_rates_; // the viscosity of the fluid mixture in the segments // TODO: it is only used to calculate the Reynolds number as we know // maybe it is not better just to store the Reynolds number? ADB segment_viscosities_; std::vector top_well_segments_; // segment volume by dt (time step) // to handle the volume effects of the segment Vector segvdt_; // technically, they are only useful for standard wells // since at the moment, we are handling both the standard // wells and the multi-segment wells under the MultisegmentWells // we need them to remove the dependency on StandardWells Vector well_perforation_densities_; Vector well_perforation_pressure_diffs_; }; } // namespace Opm #include "MultisegmentWells_impl.hpp" #endif // OPM_MULTISEGMENTWELLS_HEADER_INCLUDED