2016-04-20 09:55:17 +02:00
|
|
|
/*
|
|
|
|
|
Copyright 2016 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 <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
|
|
|
|
#include <opm/parser/eclipse/EclipseState/IOConfig/IOConfig.hpp>
|
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/ScheduleEnums.hpp>
|
|
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
2016-04-26 12:26:38 +02:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Group.hpp>
|
2016-04-20 09:55:17 +02:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
|
2016-04-26 12:26:38 +02:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/WellSet.hpp>
|
2016-04-20 09:55:17 +02:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/WellProductionProperties.hpp>
|
|
|
|
|
#include <opm/parser/eclipse/EclipseState/SummaryConfig/SummaryConfig.hpp>
|
2016-08-23 14:09:44 +02:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Grid/GridProperty.hpp>
|
2016-04-20 09:55:17 +02:00
|
|
|
#include <opm/parser/eclipse/Units/UnitSystem.hpp>
|
|
|
|
|
|
2016-10-11 14:39:15 +02:00
|
|
|
#include <opm/output/eclipse/Summary.hpp>
|
|
|
|
|
#include <opm/output/eclipse/RegionCache.hpp>
|
|
|
|
|
|
2016-04-20 09:55:17 +02:00
|
|
|
/*
|
2016-06-15 10:45:07 +02:00
|
|
|
* This class takes simulator state and parser-provided information and
|
|
|
|
|
* orchestrates ert to write simulation results as requested by the SUMMARY
|
|
|
|
|
* section in eclipse. The implementation is somewhat compact as a lot of the
|
|
|
|
|
* requested output may be similar-but-not-quite. Through various techniques
|
|
|
|
|
* the compiler writes a lot of this code for us.
|
2016-04-20 09:55:17 +02:00
|
|
|
*/
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
namespace Opm {
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
namespace {
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
using rt = data::Rates::opt;
|
|
|
|
|
using measure = UnitSystem::measure;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
/* Some numerical value with its unit tag embedded to enable caller to apply
|
2016-07-19 13:18:06 +02:00
|
|
|
* unit conversion. This removes a lot of boilerplate. ad-hoc solution to poor
|
|
|
|
|
* unit support in general.
|
2016-06-15 10:45:07 +02:00
|
|
|
*/
|
2016-07-19 13:18:06 +02:00
|
|
|
measure div_unit( measure denom, measure div ) {
|
|
|
|
|
if( denom == measure::gas_surface_rate &&
|
|
|
|
|
div == measure::liquid_surface_rate )
|
|
|
|
|
return measure::gas_oil_ratio;
|
|
|
|
|
|
|
|
|
|
if( denom == measure::liquid_surface_rate &&
|
|
|
|
|
div == measure::gas_surface_rate )
|
|
|
|
|
return measure::oil_gas_ratio;
|
|
|
|
|
|
|
|
|
|
if( denom == measure::liquid_surface_rate &&
|
|
|
|
|
div == measure::liquid_surface_rate )
|
|
|
|
|
return measure::water_cut;
|
|
|
|
|
|
|
|
|
|
if( denom == measure::liquid_surface_rate &&
|
|
|
|
|
div == measure::time )
|
|
|
|
|
return measure::liquid_surface_volume;
|
|
|
|
|
|
|
|
|
|
if( denom == measure::gas_surface_rate &&
|
|
|
|
|
div == measure::time )
|
|
|
|
|
return measure::gas_surface_volume;
|
|
|
|
|
|
|
|
|
|
return measure::identity;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
measure mul_unit( measure lhs, measure rhs ) {
|
|
|
|
|
if( lhs == rhs ) return lhs;
|
|
|
|
|
|
|
|
|
|
if( ( lhs == measure::liquid_surface_rate && rhs == measure::time ) ||
|
|
|
|
|
( rhs == measure::liquid_surface_rate && lhs == measure::time ) )
|
|
|
|
|
return measure::liquid_surface_volume;
|
|
|
|
|
|
|
|
|
|
if( ( lhs == measure::gas_surface_rate && rhs == measure::time ) ||
|
|
|
|
|
( rhs == measure::gas_surface_rate && lhs == measure::time ) )
|
|
|
|
|
return measure::gas_surface_volume;
|
|
|
|
|
|
|
|
|
|
return lhs;
|
|
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
struct quantity {
|
|
|
|
|
double value;
|
|
|
|
|
UnitSystem::measure unit;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
quantity operator+( const quantity& rhs ) const {
|
2016-07-19 13:18:06 +02:00
|
|
|
assert( this->unit == rhs.unit );
|
2016-06-15 10:45:07 +02:00
|
|
|
return { this->value + rhs.value, this->unit };
|
|
|
|
|
}
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
quantity operator*( const quantity& rhs ) const {
|
2016-07-19 13:18:06 +02:00
|
|
|
return { this->value * rhs.value, mul_unit( this->unit, rhs.unit ) };
|
2016-06-15 10:45:07 +02:00
|
|
|
}
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
quantity operator/( const quantity& rhs ) const {
|
2016-07-19 13:18:06 +02:00
|
|
|
const auto res_unit = div_unit( this->unit, rhs.unit );
|
|
|
|
|
|
|
|
|
|
if( rhs.value == 0 ) return { 0.0, res_unit };
|
|
|
|
|
return { this->value / rhs.value, res_unit };
|
2016-06-15 10:45:07 +02:00
|
|
|
}
|
2016-10-24 12:44:56 +02:00
|
|
|
|
|
|
|
|
quantity operator/( double divisor ) const {
|
|
|
|
|
if( divisor == 0 ) return { 0.0, this->unit };
|
|
|
|
|
return { this->value / divisor , this->unit };
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
quantity& operator/=( double divisor ) {
|
|
|
|
|
if( divisor == 0 )
|
|
|
|
|
this->value = 0;
|
|
|
|
|
else
|
|
|
|
|
this->value /= divisor;
|
|
|
|
|
|
|
|
|
|
return *this;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
2016-09-23 15:57:10 +02:00
|
|
|
quantity operator-( const quantity& rhs) const {
|
|
|
|
|
return { this->value - rhs.value, this->unit };
|
|
|
|
|
}
|
2016-06-15 10:45:07 +02:00
|
|
|
};
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
/*
|
|
|
|
|
* All functions must have the same parameters, so they're gathered in a struct
|
|
|
|
|
* and functions use whatever information they care about.
|
|
|
|
|
*
|
2016-07-19 08:46:09 +02:00
|
|
|
* schedule_wells are wells from the deck, provided by opm-parser. active_index
|
|
|
|
|
* is the index of the block in question. wells is simulation data.
|
2016-06-15 10:45:07 +02:00
|
|
|
*/
|
|
|
|
|
struct fn_args {
|
2016-07-19 08:46:09 +02:00
|
|
|
const std::vector< const Well* >& schedule_wells;
|
2016-06-15 10:45:07 +02:00
|
|
|
double duration;
|
|
|
|
|
size_t timestep;
|
2016-08-23 14:09:44 +02:00
|
|
|
int num;
|
2016-06-15 10:45:07 +02:00
|
|
|
const data::Wells& wells;
|
2016-08-23 14:09:44 +02:00
|
|
|
const data::Solution& state;
|
2016-10-11 14:39:15 +02:00
|
|
|
const out::RegionCache& regionCache;
|
2016-10-25 16:25:46 +02:00
|
|
|
const EclipseGrid& grid;
|
2016-06-15 10:45:07 +02:00
|
|
|
};
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
/* Since there are several enums in opm scattered about more-or-less
|
|
|
|
|
* representing the same thing. Since functions use template parameters to
|
|
|
|
|
* expand into the actual implementations we need a static way to determine
|
|
|
|
|
* what unit to tag the return value with.
|
|
|
|
|
*/
|
|
|
|
|
template< rt > constexpr
|
|
|
|
|
measure rate_unit() { return measure::liquid_surface_rate; }
|
|
|
|
|
template< Phase::PhaseEnum > constexpr
|
|
|
|
|
measure rate_unit() { return measure::liquid_surface_rate; }
|
|
|
|
|
|
|
|
|
|
template<> constexpr
|
|
|
|
|
measure rate_unit< rt::gas >() { return measure::gas_surface_rate; }
|
|
|
|
|
template<> constexpr
|
|
|
|
|
measure rate_unit< Phase::GAS >() { return measure::gas_surface_rate; }
|
|
|
|
|
|
2016-09-23 15:57:10 +02:00
|
|
|
template<> constexpr
|
|
|
|
|
measure rate_unit< rt::solvent >() { return measure::gas_surface_rate; }
|
|
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
template< rt phase, bool injection = true >
|
2016-06-15 10:45:07 +02:00
|
|
|
inline quantity rate( const fn_args& args ) {
|
|
|
|
|
double sum = 0.0;
|
|
|
|
|
|
2016-07-19 08:46:09 +02:00
|
|
|
for( const auto* sched_well : args.schedule_wells ) {
|
|
|
|
|
const auto& name = sched_well->name();
|
2016-09-25 18:31:34 +02:00
|
|
|
if( args.wells.count( name ) == 0 ) continue;
|
2016-07-18 12:02:29 +02:00
|
|
|
const auto v = args.wells.at( name ).rates.get( phase, 0.0 );
|
|
|
|
|
if( ( v > 0 ) == injection )
|
|
|
|
|
sum += v;
|
2016-06-15 10:45:07 +02:00
|
|
|
}
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
if( !injection ) sum *= -1;
|
2016-06-15 10:45:07 +02:00
|
|
|
return { sum, rate_unit< phase >() };
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
template< rt phase, bool injection = true >
|
2016-07-15 15:20:56 +02:00
|
|
|
inline quantity crate( const fn_args& args ) {
|
2016-07-19 09:45:50 +02:00
|
|
|
const quantity zero = { 0, rate_unit< phase >() };
|
2016-10-25 16:25:46 +02:00
|
|
|
// The args.num value is the literal value which will go to the
|
|
|
|
|
// NUMS array in the eclispe SMSPEC file; the values in this array
|
|
|
|
|
// are offset 1 - whereas we need to use this index here to look
|
|
|
|
|
// up a completion with offset 0.
|
|
|
|
|
const auto global_index = args.num - 1;
|
|
|
|
|
const auto active_index = args.grid.activeIndex( global_index );
|
2016-07-19 09:45:50 +02:00
|
|
|
if( args.schedule_wells.empty() ) return zero;
|
2016-07-15 15:20:56 +02:00
|
|
|
|
2016-07-19 09:45:50 +02:00
|
|
|
const auto& name = args.schedule_wells.front()->name();
|
2016-09-25 18:31:34 +02:00
|
|
|
if( args.wells.count( name ) == 0 ) return zero;
|
2016-07-18 12:02:29 +02:00
|
|
|
|
2016-07-19 09:45:50 +02:00
|
|
|
const auto& well = args.wells.at( name );
|
2016-10-24 16:52:26 +02:00
|
|
|
const auto& completion = std::find_if( well.completions.begin(),
|
|
|
|
|
well.completions.end(),
|
|
|
|
|
[=]( const data::Completion& c ) {
|
2016-10-25 16:25:46 +02:00
|
|
|
return c.index == active_index;
|
2016-10-24 16:52:26 +02:00
|
|
|
} );
|
|
|
|
|
|
|
|
|
|
if( completion == well.completions.end() ) return zero;
|
|
|
|
|
const auto v = completion->rates.get( phase, 0.0 );
|
2016-07-19 09:45:50 +02:00
|
|
|
if( ( v > 0 ) != injection ) return zero;
|
|
|
|
|
|
|
|
|
|
if( !injection ) return { -v, rate_unit< phase >() };
|
|
|
|
|
return { v, rate_unit< phase >() };
|
2016-07-15 15:20:56 +02:00
|
|
|
}
|
|
|
|
|
|
2016-09-23 15:57:10 +02:00
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
template< rt phase > inline quantity prodrate( const fn_args& args )
|
|
|
|
|
{ return rate< phase, false >( args ); }
|
|
|
|
|
template< rt phase > inline quantity injerate( const fn_args& args )
|
|
|
|
|
{ return rate< phase, true >( args ); }
|
|
|
|
|
template< rt phase > inline quantity prodcrate( const fn_args& args )
|
|
|
|
|
{ return crate< phase, false >( args ); }
|
|
|
|
|
template< rt phase > inline quantity injecrate( const fn_args& args )
|
|
|
|
|
{ return crate< phase, true >( args ); }
|
|
|
|
|
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
inline quantity bhp( const fn_args& args ) {
|
2016-07-19 13:00:25 +02:00
|
|
|
const quantity zero = { 0, measure::pressure };
|
|
|
|
|
if( args.schedule_wells.empty() ) return zero;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-09-25 18:31:34 +02:00
|
|
|
const auto p = args.wells.find( args.schedule_wells.front()->name() );
|
|
|
|
|
if( p == args.wells.end() ) return zero;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return { p->second.bhp, measure::pressure };
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
inline quantity thp( const fn_args& args ) {
|
2016-07-19 13:00:25 +02:00
|
|
|
const quantity zero = { 0, measure::pressure };
|
|
|
|
|
if( args.schedule_wells.empty() ) return zero;
|
2016-04-28 09:46:01 +02:00
|
|
|
|
2016-09-25 18:31:34 +02:00
|
|
|
const auto p = args.wells.find( args.schedule_wells.front()->name() );
|
|
|
|
|
if( p == args.wells.end() ) return zero;
|
2016-04-28 09:46:01 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return { p->second.thp, measure::pressure };
|
2016-04-28 09:46:01 +02:00
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
template< Phase::PhaseEnum phase >
|
|
|
|
|
inline quantity production_history( const fn_args& args ) {
|
2016-06-10 09:54:44 +02:00
|
|
|
/*
|
|
|
|
|
* For well data, looking up historical rates (both for production and
|
|
|
|
|
* injection) before simulation actually starts is impossible and
|
|
|
|
|
* nonsensical. We therefore default to writing zero (which is what eclipse
|
|
|
|
|
* seems to do as well). Additionally, when an input deck is parsed,
|
|
|
|
|
* timesteps and rates are structured as such:
|
|
|
|
|
*
|
|
|
|
|
* The rates observed in timestep N is denoted at timestep N-1, i.e. at the
|
|
|
|
|
* **end** of the previous timestep. Which means that what observed at
|
|
|
|
|
* timestep 1 is denoted at timestep 0, and what happens "on" timestep 0
|
|
|
|
|
* doesn't exist and would in code give an arithmetic error. We therefore
|
|
|
|
|
* special-case timestep N == 0, and for all other timesteps look up the
|
|
|
|
|
* value *reported* at N-1 which applies to timestep N.
|
|
|
|
|
*/
|
2016-06-15 10:45:07 +02:00
|
|
|
if( args.timestep == 0 ) return { 0.0, rate_unit< phase >() };
|
2016-06-10 09:54:44 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto timestep = args.timestep - 1;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
double sum = 0.0;
|
2016-07-19 08:46:09 +02:00
|
|
|
for( const Well* sched_well : args.schedule_wells )
|
|
|
|
|
sum += sched_well->production_rate( phase, timestep );
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return { sum, rate_unit< phase >() };
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
template< Phase::PhaseEnum phase >
|
|
|
|
|
inline quantity injection_history( const fn_args& args ) {
|
|
|
|
|
if( args.timestep == 0 ) return { 0.0, rate_unit< phase >() };
|
2016-05-31 14:42:59 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto timestep = args.timestep - 1;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
double sum = 0.0;
|
2016-07-19 08:46:09 +02:00
|
|
|
for( const Well* sched_well : args.schedule_wells )
|
|
|
|
|
sum += sched_well->injection_rate( phase, timestep );
|
2016-05-02 13:04:43 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return { sum, rate_unit< phase >() };
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
/*
|
|
|
|
|
* A small DSL, really poor man's function composition, to avoid massive
|
|
|
|
|
* repetition when declaring the handlers for each individual keyword. bin_op
|
|
|
|
|
* and its aliases will apply the pair of functions F and G that all take const
|
|
|
|
|
* fn_args& and return quantity, making them composable.
|
|
|
|
|
*/
|
|
|
|
|
template< typename F, typename G, typename Op >
|
|
|
|
|
struct bin_op {
|
|
|
|
|
bin_op( F fn, G gn = {} ) : f( fn ), g( gn ) {}
|
|
|
|
|
quantity operator()( const fn_args& args ) const {
|
|
|
|
|
return Op()( f( args ), g( args ) );
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
private:
|
|
|
|
|
F f;
|
|
|
|
|
G g;
|
|
|
|
|
};
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-07-19 13:18:06 +02:00
|
|
|
inline quantity duration( const fn_args& args ) {
|
|
|
|
|
return { args.duration, measure::time };
|
|
|
|
|
}
|
2016-06-15 10:45:07 +02:00
|
|
|
|
2016-10-25 16:26:14 +02:00
|
|
|
template<rt phase , bool injection>
|
|
|
|
|
quantity region_rate( const fn_args& args ) {
|
|
|
|
|
double sum = 0;
|
|
|
|
|
const auto& well_completions = args.regionCache.completions( args.num );
|
|
|
|
|
for (const auto& pair : well_completions) {
|
|
|
|
|
double rate = args.wells.get( pair.first , pair.second , phase );
|
|
|
|
|
|
|
|
|
|
// We are asking for the production rate in an injector - or
|
|
|
|
|
// opposite. We just clamp to zero.
|
|
|
|
|
if ((rate > 0) != injection)
|
|
|
|
|
rate = 0;
|
|
|
|
|
|
|
|
|
|
sum += rate;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if( injection )
|
|
|
|
|
return { sum, rate_unit< phase >() };
|
|
|
|
|
else
|
|
|
|
|
return { -sum, rate_unit< phase >() };
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
template<rt phase>
|
|
|
|
|
quantity region_production(const fn_args& args) {
|
|
|
|
|
return region_rate<phase , false>(args);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
template<rt phase>
|
|
|
|
|
quantity region_injection(const fn_args& args) {
|
|
|
|
|
return region_rate<phase , true>(args);
|
|
|
|
|
}
|
|
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
|
|
|
|
|
quantity region_sum( const fn_args& args , const std::string& keyword , UnitSystem::measure unit) {
|
2016-10-11 14:39:15 +02:00
|
|
|
const auto& cells = args.regionCache.cells( args.num );
|
2016-10-25 16:26:14 +02:00
|
|
|
if (cells.empty())
|
2016-10-24 12:44:56 +02:00
|
|
|
return { 0.0 , unit };
|
2016-08-23 14:09:44 +02:00
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
double sum = 0;
|
2016-10-06 08:16:33 +02:00
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
const std::vector<double>& sim_value = args.state.data( keyword);
|
2016-08-23 14:09:44 +02:00
|
|
|
|
|
|
|
|
for (auto cell_index : cells)
|
2016-10-24 12:44:56 +02:00
|
|
|
sum += sim_value[cell_index];
|
2016-08-23 14:09:44 +02:00
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
return { sum , unit };
|
|
|
|
|
}
|
2016-08-23 14:09:44 +02:00
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
|
|
|
|
|
quantity rpr(const fn_args& args) {
|
|
|
|
|
quantity p = region_sum( args , "PRESSURE" ,measure::pressure );
|
|
|
|
|
const auto& cells = args.regionCache.cells( args.num );
|
|
|
|
|
if (cells.size() > 0)
|
|
|
|
|
p /= cells.size();
|
|
|
|
|
return p;
|
2016-08-23 14:09:44 +02:00
|
|
|
}
|
|
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
quantity roip(const fn_args& args) {
|
|
|
|
|
return region_sum( args , "ROIP", measure::volume );
|
|
|
|
|
}
|
|
|
|
|
|
2016-10-28 12:48:12 +02:00
|
|
|
quantity rgip(const fn_args& args) {
|
|
|
|
|
return region_sum( args , "RGIP", measure::volume );
|
|
|
|
|
}
|
2016-10-24 12:44:56 +02:00
|
|
|
|
|
|
|
|
quantity roipl(const fn_args& args) {
|
|
|
|
|
return region_sum( args , "ROIPL", measure::volume );
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
quantity roipg(const fn_args& args) {
|
|
|
|
|
return region_sum( args , "ROIPG", measure::volume );
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
template< typename F, typename G >
|
2016-07-19 13:18:06 +02:00
|
|
|
auto mul( F f, G g ) -> bin_op< F, G, std::multiplies< quantity > >
|
|
|
|
|
{ return { f, g }; }
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
template< typename F, typename G >
|
2016-07-19 13:18:06 +02:00
|
|
|
auto sum( F f, G g ) -> bin_op< F, G, std::plus< quantity > >
|
|
|
|
|
{ return { f, g }; }
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
template< typename F, typename G >
|
2016-07-19 13:18:06 +02:00
|
|
|
auto div( F f, G g ) -> bin_op< F, G, std::divides< quantity > >
|
|
|
|
|
{ return { f, g }; }
|
2016-06-15 10:45:07 +02:00
|
|
|
|
2016-09-23 15:57:10 +02:00
|
|
|
template< typename F, typename G >
|
|
|
|
|
auto sub( F f, G g ) -> bin_op< F, G, std::minus< quantity > >
|
|
|
|
|
{ return { f, g }; }
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
using ofun = std::function< quantity( const fn_args& ) >;
|
|
|
|
|
|
|
|
|
|
static const std::unordered_map< std::string, ofun > funs = {
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "WWIR", injerate< rt::wat > },
|
|
|
|
|
{ "WOIR", injerate< rt::oil > },
|
|
|
|
|
{ "WGIR", injerate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "WNIR", injerate< rt::solvent > },
|
|
|
|
|
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "WWIT", mul( injerate< rt::wat >, duration ) },
|
|
|
|
|
{ "WOIT", mul( injerate< rt::oil >, duration ) },
|
|
|
|
|
{ "WGIT", mul( injerate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "WNIT", mul( injerate< rt::solvent >, duration ) },
|
2016-07-18 12:02:29 +02:00
|
|
|
|
|
|
|
|
{ "WWPR", prodrate< rt::wat > },
|
|
|
|
|
{ "WOPR", prodrate< rt::oil > },
|
|
|
|
|
{ "WGPR", prodrate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "WNPR", prodrate< rt::solvent > },
|
|
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "WLPR", sum( prodrate< rt::wat >, prodrate< rt::oil > ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "WWPT", mul( prodrate< rt::wat >, duration ) },
|
|
|
|
|
{ "WOPT", mul( prodrate< rt::oil >, duration ) },
|
|
|
|
|
{ "WGPT", mul( prodrate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "WNPT", mul( prodrate< rt::solvent >, duration ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "WLPT", mul( sum( prodrate< rt::wat >, prodrate< rt::oil > ),
|
|
|
|
|
duration ) },
|
2016-07-18 12:02:29 +02:00
|
|
|
|
|
|
|
|
{ "WWCT", div( prodrate< rt::wat >,
|
|
|
|
|
sum( prodrate< rt::wat >, prodrate< rt::oil > ) ) },
|
|
|
|
|
{ "GWCT", div( prodrate< rt::wat >,
|
|
|
|
|
sum( prodrate< rt::wat >, prodrate< rt::oil > ) ) },
|
|
|
|
|
{ "WGOR", div( prodrate< rt::gas >, prodrate< rt::oil > ) },
|
|
|
|
|
{ "GGOR", div( prodrate< rt::gas >, prodrate< rt::oil > ) },
|
|
|
|
|
{ "WGLR", div( prodrate< rt::gas >,
|
|
|
|
|
sum( prodrate< rt::wat >, prodrate< rt::oil > ) ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
{ "WBHP", bhp },
|
|
|
|
|
{ "WTHP", thp },
|
|
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "GWIR", injerate< rt::wat > },
|
|
|
|
|
{ "GOIR", injerate< rt::oil > },
|
|
|
|
|
{ "GGIR", injerate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "GNIR", injerate< rt::solvent > },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "GWIT", mul( injerate< rt::wat >, duration ) },
|
|
|
|
|
{ "GOIT", mul( injerate< rt::oil >, duration ) },
|
|
|
|
|
{ "GGIT", mul( injerate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "GNIT", mul( injerate< rt::solvent >, duration ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "GWPR", prodrate< rt::wat > },
|
|
|
|
|
{ "GOPR", prodrate< rt::oil > },
|
|
|
|
|
{ "GGPR", prodrate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "GNPR", prodrate< rt::solvent > },
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "GLPR", sum( prodrate< rt::wat >, prodrate< rt::oil > ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "GWPT", mul( prodrate< rt::wat >, duration ) },
|
|
|
|
|
{ "GOPT", mul( prodrate< rt::oil >, duration ) },
|
|
|
|
|
{ "GGPT", mul( prodrate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "GNPT", mul( prodrate< rt::solvent >, duration ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "GLPT", mul( sum( prodrate< rt::wat >, prodrate< rt::oil > ),
|
|
|
|
|
duration ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
{ "WWPRH", production_history< Phase::WATER > },
|
|
|
|
|
{ "WOPRH", production_history< Phase::OIL > },
|
|
|
|
|
{ "WGPRH", production_history< Phase::GAS > },
|
|
|
|
|
{ "WLPRH", sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
|
|
|
|
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "WWPTH", mul( production_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "WOPTH", mul( production_history< Phase::OIL >, duration ) },
|
|
|
|
|
{ "WGPTH", mul( production_history< Phase::GAS >, duration ) },
|
|
|
|
|
{ "WLPTH", mul( sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ),
|
|
|
|
|
duration ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
{ "WWIRH", injection_history< Phase::WATER > },
|
|
|
|
|
{ "WOIRH", injection_history< Phase::OIL > },
|
|
|
|
|
{ "WGIRH", injection_history< Phase::GAS > },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "WWITH", mul( injection_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "WOITH", mul( injection_history< Phase::OIL >, duration ) },
|
|
|
|
|
{ "WGITH", mul( injection_history< Phase::GAS >, duration ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
/* From our point of view, injectors don't have water cuts and div/sum will return 0.0 */
|
|
|
|
|
{ "WWCTH", div( production_history< Phase::WATER >,
|
|
|
|
|
sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) ) },
|
2016-04-26 12:26:38 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
/* We do not support mixed injections, and gas/oil is undefined when oil is
|
|
|
|
|
* zero (i.e. pure gas injector), so always output 0 if this is an injector
|
|
|
|
|
*/
|
|
|
|
|
{ "WGORH", div( production_history< Phase::GAS >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
|
|
|
|
{ "WGLRH", div( production_history< Phase::GAS >,
|
|
|
|
|
sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) ) },
|
|
|
|
|
|
|
|
|
|
{ "GWPRH", production_history< Phase::WATER > },
|
|
|
|
|
{ "GOPRH", production_history< Phase::OIL > },
|
|
|
|
|
{ "GGPRH", production_history< Phase::GAS > },
|
|
|
|
|
{ "GLPRH", sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
|
|
|
|
{ "GWIRH", injection_history< Phase::WATER > },
|
|
|
|
|
{ "GOIRH", injection_history< Phase::OIL > },
|
|
|
|
|
{ "GGIRH", injection_history< Phase::GAS > },
|
2016-10-28 12:13:10 +02:00
|
|
|
{ "GGORH", div( production_history< Phase::GAS >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
2016-10-28 12:19:07 +02:00
|
|
|
{ "GWCTH", div( production_history< Phase::WATER >,
|
|
|
|
|
sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) ) },
|
2016-06-15 10:45:07 +02:00
|
|
|
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "GWPTH", mul( production_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "GOPTH", mul( production_history< Phase::OIL >, duration ) },
|
|
|
|
|
{ "GGPTH", mul( production_history< Phase::GAS >, duration ) },
|
|
|
|
|
{ "GLPTH", mul( sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ),
|
|
|
|
|
duration ) },
|
|
|
|
|
{ "GWITH", mul( injection_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "GGITH", mul( injection_history< Phase::GAS >, duration ) },
|
2016-07-15 15:20:56 +02:00
|
|
|
|
2016-07-18 12:02:29 +02:00
|
|
|
{ "CWIR", injecrate< rt::wat > },
|
|
|
|
|
{ "CGIR", injecrate< rt::gas > },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "CWIT", mul( injecrate< rt::wat >, duration ) },
|
|
|
|
|
{ "CGIT", mul( injecrate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "CNIT", mul( injecrate< rt::solvent >, duration ) },
|
2016-07-18 12:02:29 +02:00
|
|
|
|
|
|
|
|
{ "CWPR", prodcrate< rt::wat > },
|
|
|
|
|
{ "COPR", prodcrate< rt::oil > },
|
|
|
|
|
{ "CGPR", prodcrate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
// Minus for injection rates and pluss for production rate
|
|
|
|
|
{ "CNFR", sub( prodcrate< rt::solvent>, injecrate<rt::solvent>) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "CWPT", mul( prodcrate< rt::wat >, duration ) },
|
|
|
|
|
{ "COPT", mul( prodcrate< rt::oil >, duration ) },
|
|
|
|
|
{ "CGPT", mul( prodcrate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "CNPT", mul( prodcrate< rt::solvent >, duration ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
|
|
|
|
|
{ "FWPR", prodrate< rt::wat > },
|
|
|
|
|
{ "FOPR", prodrate< rt::oil > },
|
|
|
|
|
{ "FGPR", prodrate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "FNPR", prodrate< rt::solvent > },
|
|
|
|
|
|
2016-07-18 12:57:06 +02:00
|
|
|
{ "FLPR", sum( prodrate< rt::wat >, prodrate< rt::oil > ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FWPT", mul( prodrate< rt::wat >, duration ) },
|
|
|
|
|
{ "FOPT", mul( prodrate< rt::oil >, duration ) },
|
|
|
|
|
{ "FGPT", mul( prodrate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "FNPT", mul( prodrate< rt::solvent >, duration ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FLPT", mul( sum( prodrate< rt::wat >, prodrate< rt::oil > ),
|
|
|
|
|
duration ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
|
|
|
|
|
{ "FWIR", injerate< rt::wat > },
|
|
|
|
|
{ "FOIR", injerate< rt::oil > },
|
|
|
|
|
{ "FGIR", injerate< rt::gas > },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "FNIR", injerate< rt::solvent > },
|
|
|
|
|
|
2016-07-18 12:57:06 +02:00
|
|
|
{ "FLIR", sum( injerate< rt::wat >, injerate< rt::oil > ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FWIT", mul( injerate< rt::wat >, duration ) },
|
|
|
|
|
{ "FOIT", mul( injerate< rt::oil >, duration ) },
|
|
|
|
|
{ "FGIT", mul( injerate< rt::gas >, duration ) },
|
2016-09-23 15:57:10 +02:00
|
|
|
{ "FNIT", mul( injerate< rt::solvent >, duration ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FLIT", mul( sum( injerate< rt::wat >, injerate< rt::oil > ),
|
|
|
|
|
duration ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
|
|
|
|
|
{ "FWPRH", production_history< Phase::WATER > },
|
|
|
|
|
{ "FOPRH", production_history< Phase::OIL > },
|
|
|
|
|
{ "FGPRH", production_history< Phase::GAS > },
|
|
|
|
|
{ "FLPRH", sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FWPTH", mul( production_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "FOPTH", mul( production_history< Phase::OIL >, duration ) },
|
|
|
|
|
{ "FGPTH", mul( production_history< Phase::GAS >, duration ) },
|
|
|
|
|
{ "FLPTH", mul( sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ),
|
|
|
|
|
duration ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
|
|
|
|
|
{ "FWIRH", injection_history< Phase::WATER > },
|
|
|
|
|
{ "FOIRH", injection_history< Phase::OIL > },
|
|
|
|
|
{ "FGIRH", injection_history< Phase::GAS > },
|
2016-07-19 13:18:06 +02:00
|
|
|
{ "FWITH", mul( injection_history< Phase::WATER >, duration ) },
|
|
|
|
|
{ "FOITH", mul( injection_history< Phase::OIL >, duration ) },
|
|
|
|
|
{ "FGITH", mul( injection_history< Phase::GAS >, duration ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
|
|
|
|
|
{ "FWCT", div( prodrate< rt::wat >,
|
|
|
|
|
sum( prodrate< rt::wat >, prodrate< rt::oil > ) ) },
|
|
|
|
|
{ "FGOR", div( prodrate< rt::gas >, prodrate< rt::oil > ) },
|
|
|
|
|
{ "FGLR", div( prodrate< rt::gas >,
|
|
|
|
|
sum( prodrate< rt::wat >, prodrate< rt::oil > ) ) },
|
2016-07-19 14:38:40 +02:00
|
|
|
{ "FWCTH", div( production_history< Phase::WATER >,
|
|
|
|
|
sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) ) },
|
2016-07-18 12:57:06 +02:00
|
|
|
{ "FGORH", div( production_history< Phase::GAS >,
|
|
|
|
|
production_history< Phase::OIL > ) },
|
|
|
|
|
{ "FGLRH", div( production_history< Phase::GAS >,
|
|
|
|
|
sum( production_history< Phase::WATER >,
|
|
|
|
|
production_history< Phase::OIL > ) ) },
|
2016-08-23 14:09:44 +02:00
|
|
|
|
|
|
|
|
/* Region properties */
|
|
|
|
|
{ "RPR" , rpr},
|
2016-10-24 12:44:56 +02:00
|
|
|
{ "ROIP" , roip},
|
|
|
|
|
{ "ROIPL" , roipl},
|
2016-10-25 16:26:14 +02:00
|
|
|
{ "ROIPG" , roipg},
|
2016-10-28 12:48:12 +02:00
|
|
|
{ "RGIP" , rgip},
|
2016-10-25 16:26:14 +02:00
|
|
|
{ "ROPR" , region_production<rt::oil>},
|
|
|
|
|
{ "RGPR" , region_production<rt::gas>},
|
|
|
|
|
{ "RWPR" , region_production<rt::wat>},
|
|
|
|
|
{ "ROPT" , mul( region_production<rt::oil>, duration )},
|
|
|
|
|
{ "RGPT" , mul( region_production<rt::gas>, duration )},
|
|
|
|
|
{ "RWPT" , mul( region_production<rt::wat>, duration )},
|
2016-06-15 10:45:07 +02:00
|
|
|
};
|
2016-04-26 12:26:38 +02:00
|
|
|
|
2016-10-24 12:44:56 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
inline std::vector< const Well* > find_wells( const Schedule& schedule,
|
|
|
|
|
const smspec_node_type* node,
|
|
|
|
|
size_t timestep ) {
|
2016-05-02 13:04:43 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto* name = smspec_node_get_wgname( node );
|
|
|
|
|
const auto type = smspec_node_get_var_type( node );
|
2016-04-26 12:26:38 +02:00
|
|
|
|
2016-07-15 15:20:56 +02:00
|
|
|
if( type == ECL_SMSPEC_WELL_VAR || type == ECL_SMSPEC_COMPLETION_VAR ) {
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto* well = schedule.getWell( name );
|
|
|
|
|
if( !well ) return {};
|
|
|
|
|
return { well };
|
2016-04-26 12:26:38 +02:00
|
|
|
}
|
2016-04-27 09:57:39 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
if( type == ECL_SMSPEC_GROUP_VAR ) {
|
2016-10-05 14:52:25 +02:00
|
|
|
if( !schedule.hasGroup( name ) ) return {};
|
|
|
|
|
|
|
|
|
|
const auto& group = schedule.getGroup( name );
|
2016-05-02 13:04:43 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
std::vector< const Well* > wells;
|
2016-10-05 14:52:25 +02:00
|
|
|
for( const auto& pair : group.getWells( timestep ) )
|
2016-06-15 10:45:07 +02:00
|
|
|
wells.push_back( pair.second );
|
2016-04-26 12:26:38 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return wells;
|
2016-04-26 12:26:38 +02:00
|
|
|
}
|
2016-04-27 09:57:39 +02:00
|
|
|
|
2016-07-18 12:57:06 +02:00
|
|
|
if( type == ECL_SMSPEC_FIELD_VAR )
|
|
|
|
|
return schedule.getWells();
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
return {};
|
2016-04-26 12:26:38 +02:00
|
|
|
}
|
|
|
|
|
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
namespace out {
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
class Summary::keyword_handlers {
|
|
|
|
|
public:
|
|
|
|
|
using fn = ofun;
|
|
|
|
|
std::vector< std::pair< smspec_node_type*, fn > > handlers;
|
|
|
|
|
};
|
|
|
|
|
|
2016-04-20 09:55:17 +02:00
|
|
|
Summary::Summary( const EclipseState& st, const SummaryConfig& sum ) :
|
2016-10-13 14:15:58 +02:00
|
|
|
Summary( st, sum, st.getIOConfig().fullBasePath().c_str() )
|
2016-04-20 09:55:17 +02:00
|
|
|
{}
|
|
|
|
|
|
|
|
|
|
Summary::Summary( const EclipseState& st,
|
2016-10-24 12:44:56 +02:00
|
|
|
const SummaryConfig& sum,
|
|
|
|
|
const std::string& basename ) :
|
2016-04-20 09:55:17 +02:00
|
|
|
Summary( st, sum, basename.c_str() )
|
|
|
|
|
{}
|
|
|
|
|
|
|
|
|
|
Summary::Summary( const EclipseState& st,
|
|
|
|
|
const SummaryConfig& sum,
|
|
|
|
|
const char* basename ) :
|
2016-06-15 10:45:07 +02:00
|
|
|
ecl_sum(
|
|
|
|
|
ecl_sum_alloc_writer(
|
2016-04-20 09:55:17 +02:00
|
|
|
basename,
|
2016-10-13 14:15:58 +02:00
|
|
|
st.getIOConfig().getFMTOUT(),
|
|
|
|
|
st.getIOConfig().getUNIFOUT(),
|
2016-04-20 09:55:17 +02:00
|
|
|
":",
|
2016-10-13 14:15:58 +02:00
|
|
|
st.getSchedule().posixStartTime(),
|
2016-04-20 09:55:17 +02:00
|
|
|
true,
|
2016-10-13 14:15:58 +02:00
|
|
|
st.getInputGrid().getNX(),
|
|
|
|
|
st.getInputGrid().getNY(),
|
|
|
|
|
st.getInputGrid().getNZ()
|
2016-04-20 09:55:17 +02:00
|
|
|
)
|
2016-06-15 10:45:07 +02:00
|
|
|
),
|
|
|
|
|
handlers( new keyword_handlers() )
|
2016-04-20 09:55:17 +02:00
|
|
|
{
|
2016-06-15 10:45:07 +02:00
|
|
|
/* register all keywords handlers and pair with the newly-registered ert
|
|
|
|
|
* entry.
|
|
|
|
|
*/
|
2016-04-20 09:55:17 +02:00
|
|
|
for( const auto& node : sum ) {
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto* keyword = node.keyword();
|
|
|
|
|
if( funs.find( keyword ) == funs.end() ) continue;
|
2016-04-27 10:33:51 +02:00
|
|
|
|
2016-07-19 13:00:25 +02:00
|
|
|
/* get unit strings by calling each function with dummy input */
|
|
|
|
|
const auto handle = funs.find( keyword )->second;
|
2016-08-23 14:09:44 +02:00
|
|
|
const std::vector< const Well* > dummy_wells;
|
2016-10-25 16:26:14 +02:00
|
|
|
EclipseGrid dummy_grid(1,1,1);
|
2016-10-11 14:39:15 +02:00
|
|
|
const fn_args no_args{ dummy_wells, 0, 0, 0, {} , {}, {} , dummy_grid };
|
2016-07-19 13:00:25 +02:00
|
|
|
const auto val = handle( no_args );
|
|
|
|
|
const auto* unit = st.getUnits().name( val.unit );
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
auto* nodeptr = ecl_sum_add_var( this->ecl_sum.get(), keyword,
|
2016-07-19 13:00:25 +02:00
|
|
|
node.wgname(), node.num(), unit, 0 );
|
|
|
|
|
this->handlers->handlers.emplace_back( nodeptr, handle );
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void Summary::add_timestep( int report_step,
|
2016-05-11 09:52:21 +02:00
|
|
|
double secs_elapsed,
|
2016-09-07 18:08:05 +02:00
|
|
|
const EclipseGrid& grid,
|
2016-04-20 09:55:17 +02:00
|
|
|
const EclipseState& es,
|
2016-10-11 14:39:15 +02:00
|
|
|
const RegionCache& regionCache,
|
2016-08-23 14:09:44 +02:00
|
|
|
const data::Wells& wells ,
|
|
|
|
|
const data::Solution& state) {
|
2016-05-11 09:52:21 +02:00
|
|
|
|
|
|
|
|
auto* tstep = ecl_sum_add_tstep( this->ecl_sum.get(), report_step, secs_elapsed );
|
2016-05-31 14:42:59 +02:00
|
|
|
const double duration = secs_elapsed - this->prev_time_elapsed;
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
const size_t timestep = report_step;
|
2016-10-13 14:15:58 +02:00
|
|
|
const auto& schedule = es.getSchedule();
|
2016-06-15 10:45:07 +02:00
|
|
|
|
|
|
|
|
for( auto& f : this->handlers->handlers ) {
|
2016-08-23 14:09:44 +02:00
|
|
|
const int num = smspec_node_get_num( f.first );
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto* genkey = smspec_node_get_gen_key1( f.first );
|
|
|
|
|
|
2016-07-19 08:46:09 +02:00
|
|
|
const auto schedule_wells = find_wells( schedule, f.first, timestep );
|
2016-10-11 14:39:15 +02:00
|
|
|
const auto val = f.second( { schedule_wells, duration, timestep, num, wells , state , regionCache , grid} );
|
2016-04-20 09:55:17 +02:00
|
|
|
|
2016-09-23 15:57:10 +02:00
|
|
|
const auto unit_applied_val = es.getUnits().from_si( val.unit, val.value );
|
2016-06-15 10:45:07 +02:00
|
|
|
const auto res = smspec_node_is_total( f.first ) && prev_tstep
|
|
|
|
|
? ecl_sum_tstep_get_from_key( prev_tstep, genkey ) + unit_applied_val
|
|
|
|
|
: unit_applied_val;
|
|
|
|
|
|
|
|
|
|
ecl_sum_tstep_set_from_node( tstep, f.first, res );
|
2016-04-26 12:26:38 +02:00
|
|
|
}
|
|
|
|
|
|
2016-04-20 09:55:17 +02:00
|
|
|
this->prev_tstep = tstep;
|
2016-05-31 14:42:59 +02:00
|
|
|
this->prev_time_elapsed = secs_elapsed;
|
2016-04-20 09:55:17 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
void Summary::write() {
|
|
|
|
|
ecl_sum_fwrite( this->ecl_sum.get() );
|
|
|
|
|
}
|
|
|
|
|
|
2016-06-15 10:45:07 +02:00
|
|
|
Summary::~Summary() {}
|
2016-04-20 09:55:17 +02:00
|
|
|
|
|
|
|
|
}
|
2016-06-15 10:45:07 +02:00
|
|
|
}
|