diff --git a/opm/core/wells/DynamicListEconLimited.hpp b/opm/core/wells/DynamicListEconLimited.hpp new file mode 100644 index 000000000..22123e1cb --- /dev/null +++ b/opm/core/wells/DynamicListEconLimited.hpp @@ -0,0 +1,57 @@ +/* + Copyright 2016 SINTEF ICT, Applied Mathematics. + + 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_DYNAMICLISTECONLIMITED_HPP +#define OPM_DYNAMICLISTECONLIMITED_HPP + +#include +#include + +#include + +namespace Opm +{ + + /// to handle the wells and connections voilating economic limits. + class DynamicListEconLimited + { + public: + bool anyWellEconLimited() const { + return !(m_shut_wells.empty()); + }; + + bool wellEconLimited(const std::string& well_name) const { + return std::find(m_shut_wells.begin(), m_shut_wells.end(), well_name) != m_shut_wells.end(); + }; + + void addShuttedWell(const std::string& well_name) { + // the well should not be in the list + // TODO: not sure wheter a shutted well can + // still be running through some other mechanism. + assert( !wellEconLimited(well_name) ); + + m_shut_wells.push_back(well_name); + }; + + private: + std::vector m_shut_wells; + }; + +} // namespace Opm +#endif /* OPM_DYNAMICLISTECONLIMITED_HPP */ + diff --git a/opm/core/wells/WellsManager.cpp b/opm/core/wells/WellsManager.cpp index 124d35ff1..d14035abf 100644 --- a/opm/core/wells/WellsManager.cpp +++ b/opm/core/wells/WellsManager.cpp @@ -332,12 +332,15 @@ namespace Opm const double* permeability) : w_(0), is_parallel_run_(false) { - std::vector dummy_well_potentials; + std::vector dummy_well_potentials;; + // TODO: not sure about the usage of this WellsManager constructor + // TODO: not sure whether this is the correct thing to do here. + DynamicListEconLimited dummy_list_econ_limited; init(eclipseState, timeStep, UgGridHelpers::numCells(grid), UgGridHelpers::globalCell(grid), UgGridHelpers::cartDims(grid), UgGridHelpers::dimensions(grid), UgGridHelpers::cell2Faces(grid), UgGridHelpers::beginFaceCentroids(grid), - permeability, dummy_well_potentials); + permeability, dummy_list_econ_limited, dummy_well_potentials); } diff --git a/opm/core/wells/WellsManager.hpp b/opm/core/wells/WellsManager.hpp index fc41c2fd3..c43e9d095 100644 --- a/opm/core/wells/WellsManager.hpp +++ b/opm/core/wells/WellsManager.hpp @@ -25,6 +25,7 @@ #include #include +#include #include #include @@ -86,6 +87,7 @@ namespace Opm const F2C& f2c, FC begin_face_centroids, const double* permeability, + const DynamicListEconLimited& list_econ_limited, bool is_parallel_run=false, const std::vector& well_potentials={}); @@ -155,6 +157,7 @@ namespace Opm const C2F& cell_to_faces, FC begin_face_centroids, const double* permeability, + const DynamicListEconLimited& list_econ_limited, const std::vector& well_potentials); // Disable copying and assignment. WellsManager(const WellsManager& other); @@ -178,7 +181,8 @@ namespace Opm const std::map& cartesian_to_compressed, const double* permeability, const NTG& ntg, - std::vector& wells_on_proc); + std::vector& wells_on_proc, + const DynamicListEconLimited& list_econ_limited); void addChildGroups(GroupTreeNodeConstPtr parentNode, std::shared_ptr< const Schedule > schedule, size_t timeStep, const PhaseUsage& phaseUsage); void setupGuideRates(std::vector& wells, const size_t timeStep, std::vector& well_data, std::map& well_names_to_index, diff --git a/opm/core/wells/WellsManager_impl.hpp b/opm/core/wells/WellsManager_impl.hpp index 9cfc3c92e..1ca0aa286 100644 --- a/opm/core/wells/WellsManager_impl.hpp +++ b/opm/core/wells/WellsManager_impl.hpp @@ -120,7 +120,8 @@ void WellsManager::createWellsFromSpecs(std::vector& wells, size_t const std::map& cartesian_to_compressed, const double* permeability, const NTG& ntg, - std::vector& wells_on_proc) + std::vector& wells_on_proc, + const DynamicListEconLimited& list_econ_limited) { if (dimensions != 3) { OPM_THROW(std::domain_error, @@ -144,6 +145,10 @@ void WellsManager::createWellsFromSpecs(std::vector& wells, size_t continue; } + if (list_econ_limited.wellEconLimited(well->name())) { + continue; + } + { // COMPDAT handling auto completionSet = well->getCompletions(timeStep); // shut completions and open ones stored in this process will have 1 others 0. @@ -327,13 +332,14 @@ WellsManager(const Opm::EclipseStateConstPtr eclipseState, const C2F& cell_to_faces, FC begin_face_centroids, const double* permeability, + const DynamicListEconLimited& list_econ_limited, bool is_parallel_run, const std::vector& well_potentials) : w_(0), is_parallel_run_(is_parallel_run) { init(eclipseState, timeStep, number_of_cells, global_cell, cart_dims, dimensions, - cell_to_faces, begin_face_centroids, permeability, well_potentials); + cell_to_faces, begin_face_centroids, permeability, list_econ_limited, well_potentials); } /// Construct wells from deck. @@ -348,6 +354,7 @@ WellsManager::init(const Opm::EclipseStateConstPtr eclipseState, const C2F& cell_to_faces, FC begin_face_centroids, const double* permeability, + const DynamicListEconLimited& list_econ_limited, const std::vector& well_potentials) { if (dimensions != 3) { @@ -410,7 +417,7 @@ WellsManager::init(const Opm::EclipseStateConstPtr eclipseState, dz, well_names, well_data, well_names_to_index, pu, cartesian_to_compressed, permeability, ntg, - wells_on_proc); + wells_on_proc, list_econ_limited); setupWellControls(wells, timeStep, well_names, pu, wells_on_proc);