Files
opm-common/opm/output/eclipse/AggregateWellData.hpp
Bård Skaflestad 1c21ec9d5c Restart *WEL Arrays: Add Facility for Aggregating Vectors
This commit adds a new class that aggregates the *WEL vectors for
restart purposes.  Using the data sources

    - Simulation run's "Schedule" object
    - Simulation run's "UnitSystem" object
    - Current simulation step ("sim_step")
    - WellRates object from simulator at the end of sim_step
    - SummaryState object from Summary class at the end of sim_step

this class will fill in the (presently) known elements of the IWEL
(integers), SWEL (Real/float), XWEL (double precision), and ZWEL
arrays for output to a restart file.  We distinguish contributions
of data sources extracted directly from the simulation deck from
those computed dynamically by the simulator--mostly to separate
concerns and to enable independent testing.  If this introduces too
much computational overhead in actual simulation runs we can fuse
the member functions

    AggregateWellData::captureDeclaredWellData()
    AggregateWellData::captureDynamicWellData()

to reduce the number of loops over active wells.

The overall structure of this facility is that we have a single
templated function, wellLoop(), that calls user-defined "operations"
on each active well.  Each of these operations in turn define a
subset of one of the {I,S,X,Z}WEL arrays pertaining to the
particular well.  We furthermore put implementation functions into
namespaces according to the pertinent arrays.

Many thanks to my Equinor collaborators.
2018-07-12 10:44:27 +02:00

100 lines
3.0 KiB
C++

/*
Copyright (c) 2018 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_AGGREGATE_WELL_DATA_HPP
#define OPM_AGGREGATE_WELL_DATA_HPP
#include <opm/output/eclipse/CharArrayNullTerm.hpp>
#include <opm/output/eclipse/WindowedArray.hpp>
#include <cstddef>
#include <string>
#include <vector>
namespace Opm {
class Schedule;
class SummaryState;
class UnitSystem;
} // Opm
namespace Opm { namespace data {
class WellRates;
}} // Opm::data
namespace Opm { namespace RestartIO { namespace Helpers {
class AggregateWellData
{
public:
explicit AggregateWellData(const std::vector<int>& inteHead);
void captureDeclaredWellData(const Opm::Schedule& sched,
const Opm::UnitSystem& units,
const std::size_t sim_step);
void captureDynamicWellData(const Opm::Schedule& sched,
const std::size_t sim_step,
const Opm::data::WellRates& xw,
const Opm::SummaryState& smry);
/// Retrieve Integer Well Data Array.
const std::vector<int>& getIWell() const
{
return this->iWell_.data();
}
/// Retrieve Floating-Point (Real) Well Data Array.
const std::vector<float>& getSWell() const
{
return this->sWell_.data();
}
/// Retrieve Floating-Point (Double Precision) Well Data Array.
const std::vector<double>& getXWell() const
{
return this->xWell_.data();
}
/// Retrieve Character Well Data Array.
const std::vector<CharArrayNullTerm<8>>& getZWell() const
{
return this->zWell_.data();
}
private:
/// Aggregate 'IWEL' array (Integer) for all wells.
WindowedArray<int> iWell_;
/// Aggregate 'SWEL' array (Real) for all wells.
WindowedArray<float> sWell_;
/// Aggregate 'XWEL' array (Double Precision) for all wells.
WindowedArray<double> xWell_;
/// Aggregate 'ZWEL' array (Character) for all wells.
WindowedArray<CharArrayNullTerm<8>> zWell_;
/// Maximum number of groups in model.
int nWGMax_;
};
}}} // Opm::RestartIO::Helpers
#endif // OPM_AGGREGATE_WELL_DATA_HPP