mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Significant modification of well data structures and related functions.
The following changes are made: - The SurfaceComponent enum has been removed. - Added new member Wells::number_of_phases. - The Wells::zfrac member has been replaced with comp_frac. The old zfrac always had 3 components per well (accessed according to the canonical ordering given by SurfaceComponent), the new one has number_of_phases components per well. - Changed add_well() accordingly to accept comp_frac. - Added new member WellControls::distr, giving distributions for rate controls. - All functions dealing with well controls now take Wells* and a well index instead of directly taking WellControls*. - Now append_well_controls() also takes a rate distribution argument. - Added new public function set_current_control().
This commit is contained in:
parent
743085bd16
commit
597a2cc7af
@ -35,26 +35,29 @@ extern "C" {
|
|||||||
|
|
||||||
/** Well type indicates desired/expected well behaviour. */
|
/** Well type indicates desired/expected well behaviour. */
|
||||||
enum WellType { INJECTOR, PRODUCER };
|
enum WellType { INJECTOR, PRODUCER };
|
||||||
|
|
||||||
/** Type of well control equation or inequality constraint.
|
/** Type of well control equation or inequality constraint.
|
||||||
* BHP -> Well constrained by bottom-hole pressure target.
|
* BHP -> Well constrained by bottom-hole pressure target.
|
||||||
* RATE -> Well constrained by total reservoir volume flow rate.
|
* RESERVOIR_RATE -> Well constrained by reservoir volume flow rates.
|
||||||
|
* SURFACE_RATE -> Well constrained by surface volume flow rates.
|
||||||
*/
|
*/
|
||||||
enum WellControlType { BHP , RATE };
|
enum WellControlType { BHP, RESERVOIR_RATE, SURFACE_RATE };
|
||||||
/** Canonical component names and ordering. */
|
|
||||||
enum SurfaceComponent { WATER = 0, OIL = 1, GAS = 2 };
|
|
||||||
|
|
||||||
|
|
||||||
/** Controls for a single well.
|
/** Controls for a single well.
|
||||||
|
|
||||||
* Each control specifies a well rate or bottom-hole pressure. Only
|
* Each control specifies a well rate or bottom-hole pressure. Only
|
||||||
* one control can be active at a time, indicated by current. The
|
* one control can be active at a time, indicated by current. The
|
||||||
* meaning of each control's target value depends on the control
|
* meaning of each control's target value depends on the control type:
|
||||||
* type, for BHP controls it is a pressure in Pascal, for RATE
|
* BHP -> target pressure in Pascal.
|
||||||
* controls it is a volumetric rate in cubic(meter)/second. The sign
|
* RESERVOIR_RATE -> target reservoir volume rate in cubic(meter)/second
|
||||||
* convention for RATE targets is as follows:
|
* 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 into reservoir, i.e. injecting.
|
||||||
* (-) Fluid flowing out of reservoir, i.e. producing.
|
* (-) Fluid flowing out of reservoir, i.e. producing.
|
||||||
* The active control as an equality constraint, whereas the
|
* 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
|
* non-active controls should be interpreted as inequality
|
||||||
* constraints (upper or lower bounds). For instance, a PRODUCER's
|
* constraints (upper or lower bounds). For instance, a PRODUCER's
|
||||||
* BHP constraint defines a minimum acceptable bottom-hole pressure
|
* BHP constraint defines a minimum acceptable bottom-hole pressure
|
||||||
@ -63,8 +66,10 @@ enum SurfaceComponent { WATER = 0, OIL = 1, GAS = 2 };
|
|||||||
struct WellControls
|
struct WellControls
|
||||||
{
|
{
|
||||||
int num; /** Number of controls. */
|
int num; /** Number of controls. */
|
||||||
enum WellControlType *type; /** Array of control types. */
|
enum WellControlType *type; /** Array of control types.*/
|
||||||
double *target; /** Array of control targets. */
|
double *target; /** Array of control targets */
|
||||||
|
double *distr; /** Array of rate control distributions,
|
||||||
|
Wells::number_of_phases numbers per well */
|
||||||
int current; /** Index of current active control. */
|
int current; /** Index of current active control. */
|
||||||
|
|
||||||
void *data; /** Internal management structure. */
|
void *data; /** Internal management structure. */
|
||||||
@ -76,10 +81,11 @@ struct WellControls
|
|||||||
struct Wells
|
struct Wells
|
||||||
{
|
{
|
||||||
int number_of_wells; /** Number of wells. */
|
int number_of_wells; /** Number of wells. */
|
||||||
|
int number_of_phases; /** Number of phases. */
|
||||||
|
|
||||||
enum WellType *type; /** Array of well types. */
|
enum WellType *type; /** Array of well types. */
|
||||||
double *depth_ref; /** Array of well bhp reference depths. */
|
double *depth_ref; /** Array of well bhp reference depths. */
|
||||||
double *zfrac; /** Component volume fractions for each well, size is (3*number_of_wells).
|
double *comp_frac; /** Component fractions for each well, size is (number_of_wells*number_of_phases).
|
||||||
* This is intended to be used for injection wells. For production wells
|
* This is intended to be used for injection wells. For production wells
|
||||||
* the component fractions will vary and cannot be specified a priori.
|
* the component fractions will vary and cannot be specified a priori.
|
||||||
*/
|
*/
|
||||||
@ -130,6 +136,8 @@ struct CompletionData
|
|||||||
* destroy_wells() to dispose of the Wells object and its allocated
|
* destroy_wells() to dispose of the Wells object and its allocated
|
||||||
* memory resources.
|
* memory resources.
|
||||||
*
|
*
|
||||||
|
* \param[in] nphases Number of active phases in simulation scenario.
|
||||||
|
*
|
||||||
* \param[in] nwells Expected number of wells in simulation scenario.
|
* \param[in] nwells Expected number of wells in simulation scenario.
|
||||||
* Pass zero if the total number of wells is unknown.
|
* Pass zero if the total number of wells is unknown.
|
||||||
*
|
*
|
||||||
@ -142,7 +150,7 @@ struct CompletionData
|
|||||||
* otherwise.
|
* otherwise.
|
||||||
*/
|
*/
|
||||||
struct Wells *
|
struct Wells *
|
||||||
create_wells(int nwells, int nperf);
|
create_wells(int nphases, int nwells, int nperf);
|
||||||
|
|
||||||
|
|
||||||
/**
|
/**
|
||||||
@ -156,7 +164,7 @@ create_wells(int nwells, int nperf);
|
|||||||
* \param[in] type Type of well.
|
* \param[in] type Type of well.
|
||||||
* \param[in] depth_ref Reference depth for well's BHP.
|
* \param[in] depth_ref Reference depth for well's BHP.
|
||||||
* \param[in] nperf Number of perforations.
|
* \param[in] nperf Number of perforations.
|
||||||
* \param[in] zfrac Injection fraction (three components) or NULL.
|
* \param[in] comp_frac Injection fraction array (size equal to W->number_of_phases) or NULL.
|
||||||
* \param[in] cells Grid cells in which well is perforated. Should
|
* \param[in] cells Grid cells in which well is perforated. Should
|
||||||
* ideally be track ordered.
|
* ideally be track ordered.
|
||||||
* \param[in] WI Well production index per perforation, or NULL.
|
* \param[in] WI Well production index per perforation, or NULL.
|
||||||
@ -169,7 +177,7 @@ int
|
|||||||
add_well(enum WellType type ,
|
add_well(enum WellType type ,
|
||||||
double depth_ref,
|
double depth_ref,
|
||||||
int nperf ,
|
int nperf ,
|
||||||
const double *zfrac ,
|
const double *comp_frac,
|
||||||
const int *cells ,
|
const int *cells ,
|
||||||
const double *WI ,
|
const double *WI ,
|
||||||
struct Wells *W );
|
struct Wells *W );
|
||||||
@ -181,23 +189,35 @@ add_well(enum WellType type ,
|
|||||||
* Increments ctrl->num by one if successful. Introducing a new
|
* Increments ctrl->num by one if successful. Introducing a new
|
||||||
* operational constraint does not affect the well's notion of the
|
* operational constraint does not affect the well's notion of the
|
||||||
* currently active constraint represented by ctrl->current.
|
* currently active constraint represented by ctrl->current.
|
||||||
|
* Note that *_RATE controls now require a phase distribution array
|
||||||
|
* to be associated with the control, see WellControls.
|
||||||
*
|
*
|
||||||
* \param[in] type Control type.
|
* \param[in] type Control type.
|
||||||
* \param[in] target Target value for the control.
|
* \param[in] target Target value for the control.
|
||||||
* \param[inout] ctrl Existing set of well controls.
|
* \param[in] distr Array of size W->number_of_phases or NULL.
|
||||||
|
* \param[in] well_index Index of well to receive additional control.
|
||||||
|
* \param[in,out] W Existing set of well controls.
|
||||||
* \return Non-zero (true) if successful and zero (false) otherwise.
|
* \return Non-zero (true) if successful and zero (false) otherwise.
|
||||||
*/
|
*/
|
||||||
int
|
int
|
||||||
append_well_controls(enum WellControlType type ,
|
append_well_controls(enum WellControlType type ,
|
||||||
double target,
|
double target,
|
||||||
struct WellControls *ctrl );
|
const double *distr,
|
||||||
|
int well_index,
|
||||||
|
struct Wells *W);
|
||||||
|
|
||||||
/**
|
/**
|
||||||
* Clear all controls from a well.
|
* Set the current control for a single well.
|
||||||
|
*/
|
||||||
|
void
|
||||||
|
set_current_control(int well_index, int current_control, struct Wells *W);
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Clear all controls from a single well.
|
||||||
*
|
*
|
||||||
* Does not affect the control set capacity. */
|
* Does not affect the control set capacity. */
|
||||||
void
|
void
|
||||||
clear_well_controls(struct WellControls *ctrl);
|
clear_well_controls(int well_index, struct Wells *W);
|
||||||
|
|
||||||
|
|
||||||
/**
|
/**
|
||||||
|
@ -103,6 +103,7 @@ well_controls_destroy(struct WellControls *ctrl)
|
|||||||
{
|
{
|
||||||
if (ctrl != NULL) {
|
if (ctrl != NULL) {
|
||||||
destroy_ctrl_mgmt(ctrl->data);
|
destroy_ctrl_mgmt(ctrl->data);
|
||||||
|
free (ctrl->distr);
|
||||||
free (ctrl->target);
|
free (ctrl->target);
|
||||||
free (ctrl->type);
|
free (ctrl->type);
|
||||||
}
|
}
|
||||||
@ -125,6 +126,7 @@ well_controls_create(void)
|
|||||||
ctrl->num = 0;
|
ctrl->num = 0;
|
||||||
ctrl->type = NULL;
|
ctrl->type = NULL;
|
||||||
ctrl->target = NULL;
|
ctrl->target = NULL;
|
||||||
|
ctrl->distr = NULL;
|
||||||
ctrl->current = -1;
|
ctrl->current = -1;
|
||||||
|
|
||||||
ctrl->data = create_ctrl_mgmt();
|
ctrl->data = create_ctrl_mgmt();
|
||||||
@ -141,26 +143,31 @@ well_controls_create(void)
|
|||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
static int
|
static int
|
||||||
well_controls_reserve(int nctrl, struct WellControls *ctrl)
|
well_controls_reserve(int nctrl, int nphases, struct WellControls *ctrl)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int c, ok;
|
int c, p, ok;
|
||||||
void *type, *target;
|
void *type, *target, *distr;
|
||||||
|
|
||||||
struct WellControlMgmt *m;
|
struct WellControlMgmt *m;
|
||||||
|
|
||||||
type = realloc(ctrl->type , nctrl * sizeof *ctrl->type );
|
type = realloc(ctrl->type , nctrl * sizeof *ctrl->type );
|
||||||
target = realloc(ctrl->target, nctrl * sizeof *ctrl->target);
|
target = realloc(ctrl->target, nctrl * sizeof *ctrl->target);
|
||||||
|
distr = realloc(ctrl->distr , nphases * nctrl * sizeof *ctrl->distr );
|
||||||
|
|
||||||
ok = 0;
|
ok = 0;
|
||||||
if (type != NULL) { ctrl->type = type ; ok++; }
|
if (type != NULL) { ctrl->type = type ; ok++; }
|
||||||
if (target != NULL) { ctrl->target = target; ok++; }
|
if (target != NULL) { ctrl->target = target; ok++; }
|
||||||
|
if (distr != NULL) { ctrl->distr = distr ; ok++; }
|
||||||
|
|
||||||
if (ok == 2) {
|
if (ok == 2) {
|
||||||
m = ctrl->data;
|
m = ctrl->data;
|
||||||
for (c = m->cpty; c < nctrl; c++) {
|
for (c = m->cpty; c < nctrl; c++) {
|
||||||
ctrl->type [c] = BHP;
|
ctrl->type [c] = BHP;
|
||||||
ctrl->target[c] = -1.0;
|
ctrl->target[c] = -1.0;
|
||||||
|
for (p = 0; p < nphases; ++p) {
|
||||||
|
ctrl->distr [nphases*c + p] = 0.0;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
m->cpty = nctrl;
|
m->cpty = nctrl;
|
||||||
@ -175,14 +182,16 @@ static int
|
|||||||
wells_allocate(int nwells, struct Wells *W)
|
wells_allocate(int nwells, struct Wells *W)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int ok;
|
int ok, np;
|
||||||
void *type, *depth_ref, *zfrac;
|
void *type, *depth_ref, *comp_frac;
|
||||||
void *well_connpos;
|
void *well_connpos;
|
||||||
void *ctrls;
|
void *ctrls;
|
||||||
|
|
||||||
|
np = W->number_of_phases;
|
||||||
|
|
||||||
type = realloc(W->type , 1 * nwells * sizeof *W->type);
|
type = realloc(W->type , 1 * nwells * sizeof *W->type);
|
||||||
depth_ref = realloc(W->depth_ref, 1 * nwells * sizeof *W->depth_ref);
|
depth_ref = realloc(W->depth_ref, 1 * nwells * sizeof *W->depth_ref);
|
||||||
zfrac = realloc(W->zfrac , 3 * nwells * sizeof *W->zfrac);
|
comp_frac = realloc(W->comp_frac, np * nwells * sizeof *W->comp_frac);
|
||||||
ctrls = realloc(W->ctrls , 1 * nwells * sizeof *W->ctrls);
|
ctrls = realloc(W->ctrls , 1 * nwells * sizeof *W->ctrls);
|
||||||
|
|
||||||
well_connpos = realloc(W->well_connpos,
|
well_connpos = realloc(W->well_connpos,
|
||||||
@ -191,7 +200,7 @@ wells_allocate(int nwells, struct Wells *W)
|
|||||||
ok = 0;
|
ok = 0;
|
||||||
if (type != NULL) { W->type = type ; ok++; }
|
if (type != NULL) { W->type = type ; ok++; }
|
||||||
if (depth_ref != NULL) { W->depth_ref = depth_ref ; ok++; }
|
if (depth_ref != NULL) { W->depth_ref = depth_ref ; ok++; }
|
||||||
if (zfrac != NULL) { W->zfrac = zfrac ; ok++; }
|
if (comp_frac != NULL) { W->comp_frac = comp_frac ; ok++; }
|
||||||
if (well_connpos != NULL) { W->well_connpos = well_connpos; ok++; }
|
if (well_connpos != NULL) { W->well_connpos = well_connpos; ok++; }
|
||||||
if (ctrls != NULL) { W->ctrls = ctrls ; ok++; }
|
if (ctrls != NULL) { W->ctrls = ctrls ; ok++; }
|
||||||
|
|
||||||
@ -223,7 +232,7 @@ static int
|
|||||||
initialise_new_wells(int nwells, struct Wells *W)
|
initialise_new_wells(int nwells, struct Wells *W)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int ok, w;
|
int ok, w, p;
|
||||||
|
|
||||||
struct WellMgmt *m;
|
struct WellMgmt *m;
|
||||||
|
|
||||||
@ -233,9 +242,9 @@ initialise_new_wells(int nwells, struct Wells *W)
|
|||||||
W->type [w] = PRODUCER;
|
W->type [w] = PRODUCER;
|
||||||
W->depth_ref[w] = -1.0;
|
W->depth_ref[w] = -1.0;
|
||||||
|
|
||||||
W->zfrac[3*w + WATER] = 0.0;
|
for (p = 0; p < W->number_of_phases; ++p) {
|
||||||
W->zfrac[3*w + OIL ] = 0.0;
|
W->comp_frac[W->number_of_phases*w + p] = 0.0;
|
||||||
W->zfrac[3*w + GAS ] = 0.0;
|
}
|
||||||
|
|
||||||
W->well_connpos[w + 1] = W->well_connpos[w];
|
W->well_connpos[w + 1] = W->well_connpos[w];
|
||||||
}
|
}
|
||||||
@ -322,7 +331,7 @@ wells_reserve(int nwells, int nperf, struct Wells *W)
|
|||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
struct Wells *
|
struct Wells *
|
||||||
create_wells(int nwells, int nperf)
|
create_wells(int nphases, int nwells, int nperf)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int ok;
|
int ok;
|
||||||
@ -332,10 +341,11 @@ create_wells(int nwells, int nperf)
|
|||||||
|
|
||||||
if (W != NULL) {
|
if (W != NULL) {
|
||||||
W->number_of_wells = 0;
|
W->number_of_wells = 0;
|
||||||
|
W->number_of_phases = nphases;
|
||||||
|
|
||||||
W->type = NULL;
|
W->type = NULL;
|
||||||
W->depth_ref = NULL;
|
W->depth_ref = NULL;
|
||||||
W->zfrac = NULL;
|
W->comp_frac = NULL;
|
||||||
|
|
||||||
W->well_connpos = malloc(1 * sizeof *W->well_connpos);
|
W->well_connpos = malloc(1 * sizeof *W->well_connpos);
|
||||||
W->well_cells = NULL;
|
W->well_cells = NULL;
|
||||||
@ -386,7 +396,7 @@ destroy_wells(struct Wells *W)
|
|||||||
free(W->WI);
|
free(W->WI);
|
||||||
free(W->well_cells);
|
free(W->well_cells);
|
||||||
free(W->well_connpos);
|
free(W->well_connpos);
|
||||||
free(W->zfrac);
|
free(W->comp_frac);
|
||||||
free(W->depth_ref);
|
free(W->depth_ref);
|
||||||
free(W->type);
|
free(W->type);
|
||||||
}
|
}
|
||||||
@ -417,13 +427,13 @@ int
|
|||||||
add_well(enum WellType type ,
|
add_well(enum WellType type ,
|
||||||
double depth_ref,
|
double depth_ref,
|
||||||
int nperf ,
|
int nperf ,
|
||||||
const double *zfrac , /* Injection fraction or NULL */
|
const double *comp_frac, /* Injection fraction or NULL */
|
||||||
const int *cells ,
|
const int *cells ,
|
||||||
const double *WI , /* Well index per perf (or NULL) */
|
const double *WI , /* Well index per perf (or NULL) */
|
||||||
struct Wells *W )
|
struct Wells *W )
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int ok, nw, nperf_tot, off;
|
int ok, nw, np, nperf_tot, off;
|
||||||
int nwalloc, nperfalloc;
|
int nwalloc, nperfalloc;
|
||||||
|
|
||||||
struct WellMgmt *m;
|
struct WellMgmt *m;
|
||||||
@ -459,8 +469,9 @@ add_well(enum WellType type ,
|
|||||||
W->type [nw] = type ;
|
W->type [nw] = type ;
|
||||||
W->depth_ref[nw] = depth_ref;
|
W->depth_ref[nw] = depth_ref;
|
||||||
|
|
||||||
if (zfrac != NULL) {
|
np = W->number_of_phases;
|
||||||
memcpy(W->zfrac + 3*nw, zfrac, 3 * sizeof *W->zfrac);
|
if (comp_frac != NULL) {
|
||||||
|
memcpy(W->comp_frac + np*nw, comp_frac, np * sizeof *W->comp_frac);
|
||||||
}
|
}
|
||||||
|
|
||||||
W->well_connpos[nw + 1] = off + nperf;
|
W->well_connpos[nw + 1] = off + nperf;
|
||||||
@ -473,15 +484,22 @@ add_well(enum WellType type ,
|
|||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
int
|
int
|
||||||
append_well_controls(enum WellControlType type ,
|
append_well_controls(enum WellControlType type,
|
||||||
double target,
|
double target,
|
||||||
struct WellControls *ctrl )
|
const double *distr,
|
||||||
|
int well_index,
|
||||||
|
struct Wells *W)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
int ok, alloc;
|
int ok, alloc, p, np;
|
||||||
|
struct WellControls *ctrl;
|
||||||
struct WellControlMgmt *m;
|
struct WellControlMgmt *m;
|
||||||
|
|
||||||
|
assert (W != NULL);
|
||||||
|
|
||||||
|
ctrl = W->ctrls[well_index];
|
||||||
|
np = W->number_of_phases;
|
||||||
|
|
||||||
assert (ctrl != NULL);
|
assert (ctrl != NULL);
|
||||||
|
|
||||||
m = ctrl->data;
|
m = ctrl->data;
|
||||||
@ -489,13 +507,17 @@ append_well_controls(enum WellControlType type ,
|
|||||||
|
|
||||||
if (! ok) {
|
if (! ok) {
|
||||||
alloc = alloc_size(ctrl->num, 1, m->cpty);
|
alloc = alloc_size(ctrl->num, 1, m->cpty);
|
||||||
ok = well_controls_reserve(alloc, ctrl);
|
ok = well_controls_reserve(alloc, np, ctrl);
|
||||||
}
|
}
|
||||||
|
|
||||||
if (ok) {
|
if (ok) {
|
||||||
ctrl->type [ctrl->num] = type ;
|
ctrl->type [ctrl->num] = type ;
|
||||||
ctrl->target[ctrl->num] = target;
|
ctrl->target[ctrl->num] = target;
|
||||||
|
if (distr != NULL) {
|
||||||
|
for (p = 0; p < np; ++p) {
|
||||||
|
ctrl->distr[ctrl->num * np + p] = distr[p];
|
||||||
|
}
|
||||||
|
}
|
||||||
ctrl->num += 1;
|
ctrl->num += 1;
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -505,10 +527,21 @@ append_well_controls(enum WellControlType type ,
|
|||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
void
|
void
|
||||||
clear_well_controls(struct WellControls *ctrl)
|
set_current_control(int well_index, int current_control, struct Wells *W)
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
{
|
{
|
||||||
if (ctrl != NULL) {
|
assert (W->ctrls[well_index] != NULL);
|
||||||
ctrl->num = 0;
|
|
||||||
|
W->ctrls[well_index]->current = current_control;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
/* ---------------------------------------------------------------------- */
|
||||||
|
void
|
||||||
|
clear_well_controls(int well_index, struct Wells *W)
|
||||||
|
/* ---------------------------------------------------------------------- */
|
||||||
|
{
|
||||||
|
if (W->ctrls[well_index] != NULL) {
|
||||||
|
W->ctrls[well_index]->num = 0;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
Loading…
Reference in New Issue
Block a user