mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
WellConnectionAuxiliaryModule: add a typetag independent base class
and put its code in a compile unit. allows embedding WellConnections.hpp
This commit is contained in:
@@ -23,15 +23,28 @@
|
||||
|
||||
#include <opm/models/discretization/common/baseauxiliarymodule.hh>
|
||||
|
||||
#include <opm/grid/CpGrid.hpp>
|
||||
#include <vector>
|
||||
|
||||
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
||||
namespace Dune { class CpGrid; }
|
||||
|
||||
namespace Opm
|
||||
{
|
||||
|
||||
class Schedule;
|
||||
|
||||
class WellConnectionAuxiliaryModuleGeneric
|
||||
{
|
||||
protected:
|
||||
WellConnectionAuxiliaryModuleGeneric(const Schedule& schedule,
|
||||
const Dune::CpGrid& grid);
|
||||
|
||||
std::vector<std::vector<int> > wells_;
|
||||
};
|
||||
|
||||
template<class TypeTag>
|
||||
class WellConnectionAuxiliaryModule
|
||||
: public BaseAuxiliaryModule<TypeTag>
|
||||
, private WellConnectionAuxiliaryModuleGeneric
|
||||
{
|
||||
using GlobalEqVector = GetPropType<TypeTag, Properties::GlobalEqVector>;
|
||||
using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
|
||||
@@ -43,55 +56,8 @@ public:
|
||||
|
||||
WellConnectionAuxiliaryModule(const Schedule& schedule,
|
||||
const Dune::CpGrid& grid)
|
||||
: WellConnectionAuxiliaryModuleGeneric(schedule, grid)
|
||||
{
|
||||
// Create cartesian to compressed mapping
|
||||
const auto& globalCell = grid.globalCell();
|
||||
const auto& cartesianSize = grid.logicalCartesianSize();
|
||||
|
||||
auto size = cartesianSize[0]*cartesianSize[1]*cartesianSize[2];
|
||||
|
||||
std::vector<int> cartesianToCompressed(size, -1);
|
||||
auto begin = globalCell.begin();
|
||||
|
||||
for ( auto cell = begin, end= globalCell.end(); cell != end; ++cell )
|
||||
{
|
||||
cartesianToCompressed[ *cell ] = cell - begin;
|
||||
}
|
||||
|
||||
const auto& schedule_wells = schedule.getWellsatEnd();
|
||||
wells_.reserve(schedule_wells.size());
|
||||
|
||||
// initialize the additional cell connections introduced by wells.
|
||||
for ( const auto& well : schedule_wells )
|
||||
{
|
||||
std::vector<int> compressed_well_perforations;
|
||||
// All possible completions of the well
|
||||
const auto& completionSet = well.getConnections();
|
||||
compressed_well_perforations.reserve(completionSet.size());
|
||||
|
||||
for ( size_t c=0; c < completionSet.size(); c++ )
|
||||
{
|
||||
const auto& completion = completionSet.get(c);
|
||||
int i = completion.getI();
|
||||
int j = completion.getJ();
|
||||
int k = completion.getK();
|
||||
int cart_grid_idx = i + cartesianSize[0]*(j + cartesianSize[1]*k);
|
||||
int compressed_idx = cartesianToCompressed[cart_grid_idx];
|
||||
|
||||
if ( compressed_idx >= 0 ) // Ignore completions in inactive/remote cells.
|
||||
{
|
||||
compressed_well_perforations.push_back(compressed_idx);
|
||||
}
|
||||
}
|
||||
|
||||
if ( ! compressed_well_perforations.empty() )
|
||||
{
|
||||
std::sort(compressed_well_perforations.begin(),
|
||||
compressed_well_perforations.end());
|
||||
|
||||
wells_.push_back(compressed_well_perforations);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
unsigned numDofs() const
|
||||
@@ -102,9 +68,9 @@ public:
|
||||
|
||||
void addNeighbors(std::vector<NeighborSet>& neighbors) const
|
||||
{
|
||||
for(const auto& well_perforations : wells_)
|
||||
for (const auto& well_perforations : wells_)
|
||||
{
|
||||
for(const auto& perforation : well_perforations)
|
||||
for (const auto& perforation : well_perforations)
|
||||
neighbors[perforation].insert(well_perforations.begin(),
|
||||
well_perforations.end());
|
||||
}
|
||||
@@ -119,8 +85,6 @@ public:
|
||||
}
|
||||
|
||||
private:
|
||||
|
||||
std::vector<std::vector<int> > wells_;
|
||||
};
|
||||
|
||||
} // end namespace OPM
|
||||
|
||||
Reference in New Issue
Block a user