void SimulatorBase<Implementation>::computeWellPotentials(const Wells* wells, const BlackoilState& x, const WellState& xw, std::vector<double>& well_potentials) { const int nw = wells->number_of_wells; const int np = wells->number_of_phases; well_potentials.clear(); well_potentials.resize(nw*np,0.0); for (int w = 0; w < nw; ++w) { for (int perf = wells->well_connpos[w]; perf < wells->well_connpos[w + 1]; ++perf) { const double well_cell_pressure = x.pressure()[wells->well_cells[perf]]; const double drawdown_used = well_cell_pressure - xw.perfPress()[perf]; const WellControls* ctrl = wells->ctrls[w]; const int nwc = well_controls_get_num(ctrl); //Loop over all controls until we find a BHP control //that specifies what we need... double bhp = 0.0; for (int ctrl_index=0; ctrl_index < nwc; ++ctrl_index) { if (well_controls_iget_type(ctrl, ctrl_index) == BHP) { bhp = well_controls_iget_target(ctrl, ctrl_index); } // TODO: do something for thp; } // Calculate the pressure difference in the well perforation const double dp = xw.perfPress()[perf] - xw.bhp()[w]; const double drawdown_maximum = well_cell_pressure - (bhp + dp); for (int phase = 0; phase < np; ++phase) { well_potentials[w*np + phase] += (drawdown_maximum / drawdown_used * xw.perfPhaseRates()[perf*np + phase]); } } } }
void restoreOPM_XWELKeyword(const std::string& restart_filename, int reportstep, bool unified, WellState& wellstate) { const char * keyword = "OPM_XWEL"; const char* filename = restart_filename.c_str(); ecl_file_type* file_type = ecl_file_open(filename, 0); if (file_type != NULL) { bool block_selected = unified ? ecl_file_select_rstblock_report_step(file_type , reportstep) : true; if (block_selected) { ecl_kw_type* xwel = ecl_file_iget_named_kw(file_type , keyword, 0); const double* xwel_data = ecl_kw_get_double_ptr(xwel); std::copy_n(xwel_data + wellstate.getRestartTemperatureOffset(), wellstate.temperature().size(), wellstate.temperature().begin()); std::copy_n(xwel_data + wellstate.getRestartBhpOffset(), wellstate.bhp().size(), wellstate.bhp().begin()); std::copy_n(xwel_data + wellstate.getRestartPerfPressOffset(), wellstate.perfPress().size(), wellstate.perfPress().begin()); std::copy_n(xwel_data + wellstate.getRestartPerfRatesOffset(), wellstate.perfRates().size(), wellstate.perfRates().begin()); std::copy_n(xwel_data + wellstate.getRestartWellRatesOffset(), wellstate.wellRates().size(), wellstate.wellRates().begin()); } else { std::string error_str = "Restart file " + restart_filename + " does not contain data for report step " + std::to_string(reportstep) + "!\n"; throw std::runtime_error(error_str); } ecl_file_close(file_type); } else { std::string error_str = "Restart file " + restart_filename + " not found!\n"; throw std::runtime_error(error_str); } }
/// Compute the output. void CompressibleTpfa::computeResults(BlackoilState& state, WellState& well_state) const { UnstructuredGrid* gg = const_cast<UnstructuredGrid*>(&grid_); CompletionData completion_data; completion_data.wdp = ! wellperf_wdp_.empty() ? const_cast<double*>(&wellperf_wdp_[0]) : 0; completion_data.A = ! wellperf_A_.empty() ? const_cast<double*>(&wellperf_A_[0]) : 0; completion_data.phasemob = ! wellperf_phasemob_.empty() ? const_cast<double*>(&wellperf_phasemob_[0]) : 0; cfs_tpfa_res_wells wells_tmp; wells_tmp.W = const_cast<Wells*>(wells_); wells_tmp.data = &completion_data; cfs_tpfa_res_forces forces; forces.wells = &wells_tmp; forces.src = NULL; double* wpress = ! well_state.bhp ().empty() ? & well_state.bhp ()[0] : 0; double* wflux = ! well_state.perfRates().empty() ? & well_state.perfRates()[0] : 0; cfs_tpfa_res_flux(gg, &forces, props_.numPhases(), &trans_[0], &cell_phasemob_[0], &face_phasemob_[0], &face_gravcap_[0], &state.pressure()[0], wpress, &state.faceflux()[0], wflux); cfs_tpfa_res_fpress(gg, props_.numPhases(), &htrans_[0], &face_phasemob_[0], &face_gravcap_[0], h_, &state.pressure()[0], &state.faceflux()[0], &state.facepressure()[0]); // Compute well perforation pressures (not done by the C code). if (wells_ != 0) { const int nw = wells_->number_of_wells; for (int w = 0; w < nw; ++w) { for (int j = wells_->well_connpos[w]; j < wells_->well_connpos[w+1]; ++j) { const double bhp = well_state.bhp()[w]; well_state.perfPress()[j] = bhp + wellperf_wdp_[j]; } } } }
/// Compute two-phase transport source terms from well terms. /// Note: Unlike the incompressible version of this function, /// this version computes surface volume injection rates, /// production rates are still total reservoir volumes. /// \param[in] props Fluid and rock properties. /// \param[in] wells Wells data structure. /// \param[in] well_state Well pressures and fluxes. /// \param[out] transport_src The transport source terms. They are to be interpreted depending on sign: /// (+) positive inflow of first (water) phase (surface volume), /// (-) negative total outflow of both phases (reservoir volume). void computeTransportSource(const BlackoilPropertiesInterface& props, const Wells* wells, const WellState& well_state, std::vector<double>& transport_src) { int nc = props.numCells(); transport_src.clear(); transport_src.resize(nc, 0.0); // Well contributions. if (wells) { const int nw = wells->number_of_wells; const int np = wells->number_of_phases; if (np != 2) { OPM_THROW(std::runtime_error, "computeTransportSource() requires a 2 phase case."); } std::vector<double> A(np*np); for (int w = 0; w < nw; ++w) { const double* comp_frac = wells->comp_frac + np*w; for (int perf = wells->well_connpos[w]; perf < wells->well_connpos[w + 1]; ++perf) { const int perf_cell = wells->well_cells[perf]; double perf_rate = well_state.perfRates()[perf]; if (perf_rate > 0.0) { // perf_rate is a total inflow reservoir rate, we want a surface water rate. if (wells->type[w] != INJECTOR) { std::cout << "**** Warning: crossflow in well " << w << " perf " << perf - wells->well_connpos[w] << " ignored. Reservoir rate was " << perf_rate/Opm::unit::day << " m^3/day." << std::endl; perf_rate = 0.0; } else { assert(std::fabs(comp_frac[0] + comp_frac[1] - 1.0) < 1e-6); perf_rate *= comp_frac[0]; // Water reservoir volume rate. props.matrix(1, &well_state.perfPress()[perf], comp_frac, &perf_cell, &A[0], 0); perf_rate *= A[0]; // Water surface volume rate. } } transport_src[perf_cell] += perf_rate; } } } }
void StandardWellsSolvent:: computePropertiesForWellConnectionPressures(const SolutionState& state, const WellState& xw, std::vector<double>& b_perf, std::vector<double>& rsmax_perf, std::vector<double>& rvmax_perf, std::vector<double>& surf_dens_perf) { // 1. Compute properties required by computeConnectionPressureDelta(). // Note that some of the complexity of this part is due to the function // taking std::vector<double> arguments, and not Eigen objects. const int nperf = wells().well_connpos[wells().number_of_wells]; const int nw = wells().number_of_wells; // Compute the average pressure in each well block const Vector perf_press = Eigen::Map<const V>(xw.perfPress().data(), nperf); Vector avg_press = perf_press*0; for (int w = 0; w < nw; ++w) { for (int perf = wells().well_connpos[w]; perf < wells().well_connpos[w+1]; ++perf) { const double p_above = perf == wells().well_connpos[w] ? state.bhp.value()[w] : perf_press[perf - 1]; const double p_avg = (perf_press[perf] + p_above)/2; avg_press[perf] = p_avg; } } const std::vector<int>& well_cells = wellOps().well_cells; // Use cell values for the temperature as the wells don't knows its temperature yet. const ADB perf_temp = subset(state.temperature, well_cells); // Compute b, rsmax, rvmax values for perforations. // Evaluate the properties using average well block pressures // and cell values for rs, rv, phase condition and temperature. const ADB avg_press_ad = ADB::constant(avg_press); std::vector<PhasePresence> perf_cond(nperf); for (int perf = 0; perf < nperf; ++perf) { perf_cond[perf] = (*phase_condition_)[well_cells[perf]]; } const PhaseUsage& pu = fluid_->phaseUsage(); DataBlock b(nperf, pu.num_phases); const Vector bw = fluid_->bWat(avg_press_ad, perf_temp, well_cells).value(); if (pu.phase_used[BlackoilPhases::Aqua]) { b.col(pu.phase_pos[BlackoilPhases::Aqua]) = bw; } assert((*active_)[Oil]); assert((*active_)[Gas]); const ADB perf_rv = subset(state.rv, well_cells); const ADB perf_rs = subset(state.rs, well_cells); const Vector perf_so = subset(state.saturation[pu.phase_pos[Oil]].value(), well_cells); if (pu.phase_used[BlackoilPhases::Liquid]) { const Vector bo = fluid_->bOil(avg_press_ad, perf_temp, perf_rs, perf_cond, well_cells).value(); //const V bo_eff = subset(rq_[pu.phase_pos[Oil] ].b , well_cells).value(); b.col(pu.phase_pos[BlackoilPhases::Liquid]) = bo; // const Vector rssat = fluidRsSat(avg_press, perf_so, well_cells); const Vector rssat = fluid_->rsSat(ADB::constant(avg_press), ADB::constant(perf_so), well_cells).value(); rsmax_perf.assign(rssat.data(), rssat.data() + nperf); } else { rsmax_perf.assign(0.0, nperf); } V surf_dens_copy = superset(fluid_->surfaceDensity(0, well_cells), Span(nperf, pu.num_phases, 0), nperf*pu.num_phases); for (int phase = 1; phase < pu.num_phases; ++phase) { if ( phase == pu.phase_pos[BlackoilPhases::Vapour]) { continue; // the gas surface density is added after the solvent is accounted for. } surf_dens_copy += superset(fluid_->surfaceDensity(phase, well_cells), Span(nperf, pu.num_phases, phase), nperf*pu.num_phases); } if (pu.phase_used[BlackoilPhases::Vapour]) { // Unclear wether the effective or the pure values should be used for the wells // the current usage of unmodified properties values gives best match. //V bg_eff = subset(rq_[pu.phase_pos[Gas]].b,well_cells).value(); Vector bg = fluid_->bGas(avg_press_ad, perf_temp, perf_rv, perf_cond, well_cells).value(); Vector rhog = fluid_->surfaceDensity(pu.phase_pos[BlackoilPhases::Vapour], well_cells); // to handle solvent related if (has_solvent_) { const Vector bs = solvent_props_->bSolvent(avg_press_ad,well_cells).value(); //const V bs_eff = subset(rq_[solvent_pos_].b,well_cells).value(); // number of cells const int nc = state.pressure.size(); const ADB zero = ADB::constant(Vector::Zero(nc)); const ADB& ss = state.solvent_saturation; const ADB& sg = ((*active_)[ Gas ] ? state.saturation[ pu.phase_pos[ Gas ] ] : zero); Selector<double> zero_selector(ss.value() + sg.value(), Selector<double>::Zero); Vector F_solvent = subset(zero_selector.select(ss, ss / (ss + sg)),well_cells).value(); Vector injectedSolventFraction = Eigen::Map<const Vector>(&xw.solventFraction()[0], nperf); Vector isProducer = Vector::Zero(nperf); Vector ones = Vector::Constant(nperf,1.0); for (int w = 0; w < nw; ++w) { if(wells().type[w] == PRODUCER) { for (int perf = wells().well_connpos[w]; perf < wells().well_connpos[w+1]; ++perf) { isProducer[perf] = 1; } } } F_solvent = isProducer * F_solvent + (ones - isProducer) * injectedSolventFraction; bg = bg * (ones - F_solvent); bg = bg + F_solvent * bs; const Vector& rhos = solvent_props_->solventSurfaceDensity(well_cells); rhog = ( (ones - F_solvent) * rhog ) + (F_solvent * rhos); } b.col(pu.phase_pos[BlackoilPhases::Vapour]) = bg; surf_dens_copy += superset(rhog, Span(nperf, pu.num_phases, pu.phase_pos[BlackoilPhases::Vapour]), nperf*pu.num_phases); // const Vector rvsat = fluidRvSat(avg_press, perf_so, well_cells); const Vector rvsat = fluid_->rvSat(ADB::constant(avg_press), ADB::constant(perf_so), well_cells).value(); rvmax_perf.assign(rvsat.data(), rvsat.data() + nperf); } else { rvmax_perf.assign(0.0, nperf); } // b and surf_dens_perf is row major, so can just copy data. b_perf.assign(b.data(), b.data() + nperf * pu.num_phases); surf_dens_perf.assign(surf_dens_copy.data(), surf_dens_copy.data() + nperf * pu.num_phases); }