adding the WellInteface

will be the base class for different well models.
This commit is contained in:
Kai Bao 2017-06-15 11:34:07 +02:00
parent 060631c5dc
commit 910fe0318c
3 changed files with 412 additions and 0 deletions

View File

@ -49,6 +49,7 @@ list (APPEND MAIN_SOURCE_FILES
opm/autodiff/VFPInjProperties.cpp
opm/autodiff/WellMultiSegment.cpp
opm/autodiff/MultisegmentWells.cpp
opm/autodiff/WellInterface.cpp
opm/autodiff/MissingFeatures.cpp
opm/polymer/PolymerState.cpp
opm/polymer/PolymerBlackoilState.cpp
@ -239,6 +240,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/autodiff/WellHelpers.hpp
opm/autodiff/StandardWells.hpp
opm/autodiff/StandardWells_impl.hpp
opm/autodiff/WellInterface.hpp
opm/autodiff/StandardWellsDense.hpp
opm/autodiff/StandardWellsSolvent.hpp
opm/autodiff/StandardWellsSolvent_impl.hpp

View File

@ -0,0 +1,246 @@
/*
Copyright 2017 SINTEF ICT, Applied Mathematics.
Copyright 2017 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 <http://www.gnu.org/licenses/>.
*/
#include "config.h"
#include <opm/autodiff/WellInterface.hpp>
namespace Opm
{
WellInterface::
WellInterface(const Well* well, const size_t time_step, const Wells* wells)
{
// TODO: trying to use wells struct as little as possible here, be prepared to
// remove the wells struct in future
const std::string& well_name = well->name();
// looking for the location of the well in the wells struct
int index_well;
for (index_well = 0; index_well < wells->number_of_wells; ++index_well) {
if (well_name == std::string(wells->name[index_well])) {
break;
}
}
// should not enter the constructor if the well does not exist in the wells struct
// here, just another assertion.
assert(index_well != wells->number_of_wells);
name_ = well_name;
index_of_well_ = index_well;
well_type_ = wells->type[index_well];
number_of_phases_ = wells->number_of_phases;
// copying the comp_frac
{
comp_frac_.resize(number_of_phases_);
const int index_begin = index_well * number_of_phases_;
std::copy(wells->comp_frac + index_begin,
wells->comp_frac + index_begin + number_of_phases_, comp_frac_.begin() );
}
well_controls_ = wells->ctrls[index_well];
// perforations related
{
const int perf_index_begin = wells->well_connpos[index_well];
const int perf_index_end = wells->well_connpos[index_well + 1];
number_of_perforations_ = perf_index_end - perf_index_begin;
well_cell_.resize(number_of_perforations_);
std::copy(wells->well_cells + perf_index_begin,
wells->well_cells + perf_index_end,
well_cell_.begin() );
well_index_.resize(number_of_perforations_);
std::copy(wells->WI + perf_index_begin,
wells->WI + perf_index_end,
well_index_.begin() );
// TODO: not sure about the processing of depth for perforations here
// Will revisit here later. There are different ways and the definition for different wells
// can be different, it is possible that we need to remove this from the WellInterface
perf_depth_.resize(number_of_perforations_, 0.);
const auto& completion_set = well->getCompletions(time_step);
for (int i = 0; i < number_of_perforations_; ++i) {
perf_depth_[i] = completion_set.get(i).getCenterDepth();
}
}
}
void
WellInterface::
init(const PhaseUsage* phase_usage_arg,
const std::vector<bool>* active_arg,
const VFPProperties* vfp_properties_arg,
const double gravity_arg)
{
phase_usage_ = phase_usage_arg;
active_ = active_arg;
vfp_properties_ = vfp_properties_arg;
gravity_ = gravity_arg;
}
const std::string&
WellInterface::
name() const
{
return name_;
}
int
WellInterface::
indexOfWell() const
{
return index_of_well_;
}
WellType
WellInterface::
wellType() const
{
return well_type_;
}
int
WellInterface::
numberOfPhases() const
{
return number_of_phases_;
}
const std::vector<double>&
WellInterface::
compFrac() const
{
return comp_frac_;
}
WellControls*
WellInterface::
wellControls() const
{
return well_controls_;
}
int
WellInterface::
numberOfPerforations() const
{
return number_of_perforations_;
}
const std::vector<double>&
WellInterface::
wellIndex() const
{
return well_index_;
}
const std::vector<double>&
WellInterface::
perfDepth() const
{
return perf_depth_;
}
const std::vector<int>&
WellInterface::
wellCells() const
{
return well_cell_;
}
const std::vector<bool>&
WellInterface::
active() const
{
assert(active_);
return *active_;
}
const PhaseUsage&
WellInterface::
phaseUsage() const
{
assert(phase_usage_);
return *phase_usage_;
}
}

View File

@ -0,0 +1,164 @@
/*
Copyright 2017 SINTEF ICT, Applied Mathematics.
Copyright 2017 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 <http://www.gnu.org/licenses/>.
*/
#ifndef OPM_WELLINTERFACE_HEADER_INCLUDED
#define OPM_WELLINTERFACE_HEADER_INCLUDED
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
#include <opm/core/wells.h>
#include <opm/core/well_controls.h>
#include <opm/core/props/BlackoilPhases.hpp>
#include <opm/autodiff/VFPProperties.hpp>
#include <opm/autodiff/VFPInjProperties.hpp>
#include <opm/autodiff/VFPProdProperties.hpp>
#include <opm/autodiff/WellHelpers.hpp>
#include <opm/autodiff/WellStateFullyImplicitBlackoilDense.hpp>
#include <string>
#include <memory>
#include <vector>
#include <cassert>
namespace Opm
{
class SimulatorFullyImplicitBlackoilEbos;
class WellInterface
{
public:
// TODO: Simulator will probably enter template parameter through TypeTag later
using Simulator = SimulatorFullyImplicitBlackoilEbos;
using WellState = WellStateFullyImplicitBlackoilDense;
/// Constructor
WellInterface(const Well* well, const size_t time_step, const Wells* wells);
/// Well name.
const std::string& name() const;
/// The index of the well in Wells struct
// It is used to locate the inforation in Wells and also WellState for now.
int indexOfWell() const;
/// Well type, INJECTOR or PRODUCER.
WellType wellType() const;
/// number of phases
int numberOfPhases() const;
/// Component fractions for each phase for the well
const std::vector<double>& compFrac() const;
/// Well controls
// TODO: later to see whether we need to return const.
WellControls* wellControls() const;
/// Number of the perforations
int numberOfPerforations() const;
/// Well productivity index for each perforation.
const std::vector<double>& wellIndex() const;
/// Depth of perforations
const std::vector<double>& perfDepth() const;
/// Indices of the grid cells/blocks that perforations are completed within.
const std::vector<int>& wellCells() const;
// TODO: the following function should be able to be removed by refactoring the well class
// It is probably only needed for StandardWell
/// the densities of the fluid in each perforation
virtual const std::vector<double>& perfDensities() const = 0;
virtual std::vector<double>& perfDensities() = 0;
/// the pressure difference between different perforations
virtual const std::vector<double>& perfPressureDiffs() const = 0;
virtual std::vector<double>& perfPressureDiffs() = 0;
virtual void assembleWellEq(Simulator& ebos_simulator,
const double dt,
WellState& well_state,
bool only_wells) = 0;
void init(const PhaseUsage* phase_usage_arg,
const std::vector<bool>* active_arg,
const VFPProperties* vfp_properties_arg,
const double gravity_arg);
// TODO: temporary
virtual void setWellVariables(const WellState& well_state) = 0;
const std::vector<bool>& active() const;
const PhaseUsage& phaseUsage() const;
protected:
// well name
std::string name_;
// the index of well in Wells struct
int index_of_well_;
// well type
// INJECTOR or PRODUCER
enum WellType well_type_;
// number of phases
int number_of_phases_;
// component fractions for each well
// typically, it should apply to injection wells
std::vector<double> comp_frac_;
// controls for this well
// TODO: later will check whehter to let it stay with pointer
struct WellControls* well_controls_;
// number of the perforations for this well
int number_of_perforations_;
// well index for each perforation
std::vector<double> well_index_;
// depth for each perforation
std::vector<double> perf_depth_;
// cell index for each well perforation
std::vector<int> well_cell_;
const PhaseUsage* phase_usage_;
const std::vector<bool>* active_;
const VFPProperties* vfp_properties_;
double gravity_;
};
}
#endif // OPM_WELLINTERFACE_HEADER_INCLUDED