Skip to content

Commit

Permalink
Remove Unused Object
Browse files Browse the repository at this point in the history
  • Loading branch information
bska committed Dec 19, 2024
1 parent 45e9590 commit f53428c
Show file tree
Hide file tree
Showing 2 changed files with 32 additions and 26 deletions.
57 changes: 31 additions & 26 deletions opm/simulators/wells/WellInterfaceGeneric.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -20,6 +20,7 @@
*/

#include <config.h>

#include <opm/simulators/wells/WellInterfaceGeneric.hpp>

#include <opm/common/ErrorMacros.hpp>
Expand All @@ -33,7 +34,9 @@
#include <opm/input/eclipse/Schedule/Well/WellPolymerProperties.hpp>
#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
#include <opm/input/eclipse/Schedule/Well/WVFPEXP.hpp>

#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>

#include <opm/simulators/wells/PerforationData.hpp>
#include <opm/simulators/wells/ParallelWellInfo.hpp>
#include <opm/simulators/wells/VFPHelpers.hpp>
Expand All @@ -45,10 +48,13 @@

#include <fmt/format.h>

#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstddef>
#include <stdexcept>
#include <string>
#include <vector>

namespace Opm {

Expand Down Expand Up @@ -670,33 +676,32 @@ void WellInterfaceGeneric<Scalar>::resetWellOperability()
}

template<class Scalar>
void WellInterfaceGeneric<Scalar>::addPerforations(const std::vector<std::tuple<int,double,double>>& perfs){
int newperf = perfs.size();

for(auto& perf : perfs){
const auto [cell, WI, depth] = perf;
auto it = std::find(well_cells_.begin(), well_cells_.end(), cell);
if(it != well_cells_.end()){
// if perforation to cell already exist just add
size_t ind = it- well_cells_.begin();
well_index_[ind] += WI;
}else{
well_cells_.push_back(cell);
well_index_.push_back(WI);
perf_depth_.push_back(depth);
// not needed
double nan = std::nan("1");
perf_rep_radius_.push_back(nan);
perf_length_.push_back(nan);
bore_diameters_.push_back(nan);
saturation_table_number_.push_back(saturation_table_number_[0]);// for now used the saturation table for the first cell
number_of_perforations_ += 1;
}
// completions_;
// inj_multiplier_;
// prev_inj_multiplier_;
// inj_fc_multiplier_;
void WellInterfaceGeneric<Scalar>::addPerforations(const std::vector<std::tuple<int,double,double>>& perfs)
{
for (const auto& [cell, WI, depth] : perfs) {
auto it = std::find(well_cells_.begin(), well_cells_.end(), cell);
if (it != well_cells_.end()) {
// If perforation to cell already exists, just add contribution.
const auto ind = std::distance(well_cells_.begin(), it);
well_index_[ind] += WI;
}
else {
well_cells_.push_back(cell);
well_index_.push_back(WI);
perf_depth_.push_back(depth);

// Not strictly needed.
const double nan = std::nan("1");
perf_rep_radius_.push_back(nan);
perf_length_.push_back(nan);
bore_diameters_.push_back(nan);

// For now use the saturation table for the first cell.
saturation_table_number_.push_back(saturation_table_number_[0]);

number_of_perforations_ += 1;
}
}
}

template<class Scalar>
Expand Down
1 change: 1 addition & 0 deletions opm/simulators/wells/WellInterfaceGeneric.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -190,6 +190,7 @@ class WellInterfaceGeneric {
void resetWellOperability();

void addPerforations(const std::vector<std::tuple<int,double,double>>& perfs);

protected:
bool getAllowCrossFlow() const;

Expand Down

0 comments on commit f53428c

Please sign in to comment.