2014-01-05 07:47:15 -06:00
|
|
|
/*===========================================================================
|
|
|
|
//
|
|
|
|
// File: newwells.c
|
|
|
|
//
|
|
|
|
// Created: 2012-02-03 11:28:40+0100
|
|
|
|
//
|
|
|
|
// Authors: Knut-Andreas Lie <Knut-Andreas.Lie@sintef.no>
|
|
|
|
// Jostein R. Natvig <Jostein.R.Natvig@sintef.no>
|
|
|
|
// Halvor M. Nilsen <HalvorMoll.Nilsen@sintef.no>
|
|
|
|
// Atgeirr F. Rasmussen <atgeirr@sintef.no>
|
|
|
|
// Xavier Raynaud <Xavier.Raynaud@sintef.no>
|
|
|
|
// Bård Skaflestad <Bard.Skaflestad@sintef.no>
|
|
|
|
//
|
|
|
|
//==========================================================================*/
|
|
|
|
|
|
|
|
|
|
|
|
/*
|
|
|
|
Copyright 2012 SINTEF ICT, Applied Mathematics.
|
|
|
|
Copyright 2012 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/core/well_controls.h>
|
2014-01-05 09:22:30 -06:00
|
|
|
#include <assert.h>
|
2014-01-05 07:47:15 -06:00
|
|
|
#include <string.h>
|
|
|
|
#include <stdlib.h>
|
2014-01-05 09:22:30 -06:00
|
|
|
#include <stdio.h>
|
2014-01-05 07:47:15 -06:00
|
|
|
|
2014-01-06 08:13:32 -06:00
|
|
|
/**
|
|
|
|
* Controls for a single well.
|
|
|
|
* Each control specifies a well rate or bottom-hole pressure. Only
|
|
|
|
* one control can be active at a time, indicated by current. The
|
|
|
|
* meaning of each control's target value depends on the control type:
|
|
|
|
*
|
|
|
|
* - BHP -> target pressure in Pascal.
|
|
|
|
* - RESERVOIR_RATE -> target reservoir volume rate in cubic(meter)/second
|
|
|
|
* - SURFACE_RATE -> target surface volume rate in cubic(meter)/second
|
|
|
|
*
|
|
|
|
* The sign convention for RATE targets is as follows:
|
|
|
|
*
|
|
|
|
* - (+) Fluid flowing into reservoir, i.e. injecting.
|
|
|
|
* - (-) Fluid flowing out of reservoir, i.e. producing.
|
|
|
|
*
|
|
|
|
* For *_RATE controls, the distribution of phases used for the control
|
|
|
|
* is also needed. For example, a total rate control should have 1.0
|
|
|
|
* for each phase, whereas a control on oil rate should have 1.0 for
|
|
|
|
* the oil phase and 0.0 for the rest. For BHP controls, this is unused.
|
|
|
|
* The active control acts as an equality constraint, whereas the
|
|
|
|
* non-active controls should be interpreted as inequality
|
|
|
|
* constraints (upper or lower bounds). For instance, a PRODUCER's
|
|
|
|
* BHP constraint defines a minimum acceptable bottom-hole pressure
|
|
|
|
* value for the well.
|
|
|
|
*/
|
|
|
|
|
|
|
|
struct WellControls
|
|
|
|
{
|
|
|
|
/**
|
|
|
|
* Number of controls.
|
|
|
|
*/
|
|
|
|
int num;
|
|
|
|
|
|
|
|
int number_of_phases;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Array of control types.
|
|
|
|
*/
|
|
|
|
enum WellControlType *type;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Array of control targets.
|
|
|
|
*/
|
|
|
|
double *target;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Array of rate control distributions,
|
|
|
|
* <CODE>number_of_phases</CODE> numbers for each control
|
|
|
|
*/
|
|
|
|
double *distr;
|
|
|
|
|
|
|
|
/**
|
|
|
|
* Index of current active control.
|
|
|
|
*/
|
|
|
|
int current;
|
|
|
|
|
|
|
|
/*
|
|
|
|
The capacity allocated.
|
|
|
|
*/
|
|
|
|
int cpty;
|
|
|
|
};
|
2014-01-05 07:47:15 -06:00
|
|
|
|
|
|
|
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
void
|
|
|
|
well_controls_destroy(struct WellControls *ctrl)
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
{
|
|
|
|
if (ctrl != NULL) {
|
|
|
|
free (ctrl->distr);
|
|
|
|
free (ctrl->target);
|
|
|
|
free (ctrl->type);
|
|
|
|
}
|
|
|
|
|
|
|
|
free(ctrl);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
struct WellControls *
|
|
|
|
well_controls_create(void)
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
{
|
|
|
|
struct WellControls *ctrl;
|
|
|
|
|
|
|
|
ctrl = malloc(1 * sizeof *ctrl);
|
|
|
|
|
|
|
|
if (ctrl != NULL) {
|
|
|
|
/* Initialise empty control set */
|
|
|
|
ctrl->num = 0;
|
|
|
|
ctrl->number_of_phases = 0;
|
|
|
|
ctrl->type = NULL;
|
|
|
|
ctrl->target = NULL;
|
|
|
|
ctrl->distr = NULL;
|
|
|
|
ctrl->current = -1;
|
2014-01-05 07:54:56 -06:00
|
|
|
ctrl->cpty = 0;
|
2014-01-05 07:47:15 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
return ctrl;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
/* ---------------------------------------------------------------------- */
|
2014-01-05 09:22:30 -06:00
|
|
|
static int
|
|
|
|
well_controls_reserve(int nctrl, struct WellControls *ctrl)
|
2014-01-05 07:47:15 -06:00
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
{
|
|
|
|
int c, p, ok;
|
|
|
|
void *type, *target, *distr;
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
type = realloc(ctrl->type , nctrl * 1 * sizeof *ctrl->type );
|
|
|
|
target = realloc(ctrl->target, nctrl * 1 * sizeof *ctrl->target);
|
|
|
|
distr = realloc(ctrl->distr , nctrl * ctrl->number_of_phases * sizeof *ctrl->distr );
|
2014-01-05 07:47:15 -06:00
|
|
|
|
|
|
|
ok = 0;
|
|
|
|
if (type != NULL) { ctrl->type = type ; ok++; }
|
|
|
|
if (target != NULL) { ctrl->target = target; ok++; }
|
|
|
|
if (distr != NULL) { ctrl->distr = distr ; ok++; }
|
|
|
|
|
|
|
|
if (ok == 3) {
|
2014-01-05 07:54:56 -06:00
|
|
|
for (c = ctrl->cpty; c < nctrl; c++) {
|
2014-01-05 07:47:15 -06:00
|
|
|
ctrl->type [c] = BHP;
|
|
|
|
ctrl->target[c] = -1.0;
|
|
|
|
}
|
|
|
|
|
2014-01-05 07:54:56 -06:00
|
|
|
for (p = ctrl->cpty * ctrl->number_of_phases; p < nctrl * ctrl->number_of_phases; ++p) {
|
2014-01-05 07:47:15 -06:00
|
|
|
ctrl->distr[ p ] = 0.0;
|
|
|
|
}
|
|
|
|
|
2014-01-05 07:54:56 -06:00
|
|
|
ctrl->cpty = nctrl;
|
2014-01-05 07:47:15 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
return ok == 3;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
int well_controls_get_num(const struct WellControls *ctrl) {
|
|
|
|
return ctrl->num;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
int well_controls_get_current( const struct WellControls * ctrl) {
|
|
|
|
return ctrl->current;
|
|
|
|
}
|
|
|
|
|
|
|
|
void
|
|
|
|
well_controls_set_current( struct WellControls * ctrl, int current) {
|
|
|
|
ctrl->current = current;
|
|
|
|
}
|
|
|
|
|
2014-01-06 05:06:36 -06:00
|
|
|
void
|
|
|
|
well_controls_invert_current( struct WellControls * ctrl ) {
|
|
|
|
ctrl->current = ~ctrl->current;
|
|
|
|
}
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
|
|
|
|
enum WellControlType
|
|
|
|
well_controls_iget_type(const struct WellControls * ctrl, int control_index) {
|
|
|
|
return ctrl->type[control_index];
|
|
|
|
}
|
|
|
|
|
2014-01-06 07:40:03 -06:00
|
|
|
|
|
|
|
enum WellControlType
|
|
|
|
well_controls_get_current_type(const struct WellControls * ctrl) {
|
|
|
|
return well_controls_iget_type( ctrl , ctrl->current);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-06 05:06:36 -06:00
|
|
|
void
|
|
|
|
well_controls_iset_type( struct WellControls * ctrls , int control_index , enum WellControlType type) {
|
|
|
|
ctrls->type[control_index] = type;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
double
|
|
|
|
well_controls_iget_target(const struct WellControls * ctrl, int control_index) {
|
|
|
|
return ctrl->target[control_index];
|
|
|
|
}
|
|
|
|
|
2014-01-06 07:40:03 -06:00
|
|
|
double
|
|
|
|
well_controls_get_current_target(const struct WellControls * ctrl) {
|
|
|
|
return ctrl->target[ctrl->current];
|
|
|
|
}
|
|
|
|
|
2014-01-06 05:06:36 -06:00
|
|
|
void
|
|
|
|
well_controls_iset_target(struct WellControls * ctrl, int control_index , double target) {
|
|
|
|
ctrl->target[control_index] = target;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
const double *
|
|
|
|
well_controls_iget_distr(const struct WellControls * ctrl, int control_index) {
|
|
|
|
int offset = control_index * ctrl->number_of_phases;
|
|
|
|
return &ctrl->distr[offset];
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-06 07:40:03 -06:00
|
|
|
const double *
|
|
|
|
well_controls_get_current_distr(const struct WellControls * ctrl) {
|
|
|
|
return well_controls_iget_distr( ctrl , ctrl->current );
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
2014-01-06 05:22:15 -06:00
|
|
|
void
|
|
|
|
well_controls_iset_distr(const struct WellControls * ctrl, int control_index, const double * distr) {
|
|
|
|
int offset = control_index * ctrl->number_of_phases;
|
|
|
|
for (int p=0; p < ctrl->number_of_phases; p++)
|
|
|
|
ctrl->distr[offset + p] = distr[p];
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
void
|
|
|
|
well_controls_assert_number_of_phases(struct WellControls * ctrl , int number_of_phases) {
|
|
|
|
if (ctrl->num == 0)
|
|
|
|
ctrl->number_of_phases = number_of_phases;
|
|
|
|
|
|
|
|
assert( ctrl->number_of_phases == number_of_phases );
|
|
|
|
}
|
|
|
|
|
|
|
|
void
|
|
|
|
well_controls_clear(struct WellControls * ctrl) {
|
|
|
|
ctrl->num = 0;
|
|
|
|
ctrl->number_of_phases = 0;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
int
|
|
|
|
well_controls_add_new(enum WellControlType type , double target , const double * distr , struct WellControls * ctrl) {
|
|
|
|
if (ctrl->num == ctrl->cpty) {
|
|
|
|
int new_cpty = 2*ctrl->cpty;
|
|
|
|
if (new_cpty == ctrl->num)
|
|
|
|
new_cpty += 1;
|
|
|
|
|
|
|
|
if (!well_controls_reserve( new_cpty , ctrl))
|
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|
2014-01-06 05:24:32 -06:00
|
|
|
well_controls_iset_type( ctrl , ctrl->num , type);
|
|
|
|
well_controls_iset_target( ctrl , ctrl->num , target);
|
2014-01-05 09:22:30 -06:00
|
|
|
|
2014-01-06 05:24:32 -06:00
|
|
|
if (distr != NULL)
|
|
|
|
well_controls_iset_distr( ctrl , ctrl->num , distr);
|
|
|
|
|
2014-01-05 09:22:30 -06:00
|
|
|
ctrl->num += 1;
|
|
|
|
return 1;
|
|
|
|
}
|
2014-01-05 07:47:15 -06:00
|
|
|
|
|
|
|
|
|
|
|
bool
|
|
|
|
well_controls_equal(const struct WellControls *ctrls1, const struct WellControls *ctrls2)
|
|
|
|
/* ---------------------------------------------------------------------- */
|
|
|
|
{
|
|
|
|
bool are_equal = true;
|
|
|
|
are_equal = (ctrls1->num == ctrls2->num);
|
2014-01-06 12:22:39 -06:00
|
|
|
are_equal = are_equal && (ctrls1->number_of_phases == ctrls2->number_of_phases);
|
2014-01-05 07:47:15 -06:00
|
|
|
if (!are_equal) {
|
|
|
|
return are_equal;
|
|
|
|
}
|
|
|
|
|
2014-01-06 12:22:39 -06:00
|
|
|
are_equal = are_equal && (memcmp(ctrls1->type, ctrls2->type, ctrls1->num * sizeof *ctrls1->type ) == 0);
|
|
|
|
are_equal = are_equal && (memcmp(ctrls1->target, ctrls2->target, ctrls1->num * sizeof *ctrls1->target ) == 0);
|
|
|
|
are_equal = are_equal && (memcmp(ctrls1->distr, ctrls2->distr, ctrls1->num * ctrls1->number_of_phases * sizeof *ctrls1->distr ) == 0);
|
|
|
|
are_equal = are_equal && (ctrls1->cpty == ctrls2->cpty);
|
2014-01-05 07:47:15 -06:00
|
|
|
|
|
|
|
return are_equal;
|
|
|
|
}
|
2014-01-05 09:22:30 -06:00
|
|
|
|