/// @brief Computes injected and produced volumes of all phases. /// Note 1: assumes that only the first phase is injected. /// Note 2: assumes that transport has been done with an /// implicit method, i.e. that the current state /// gives the mobilities used for the preceding timestep. /// @param[in] props fluid and rock properties. /// @param[in] s saturation values (for all P phases) /// @param[in] src if < 0: total outflow, if > 0: first phase inflow. /// @param[in] dt timestep used /// @param[out] injected must point to a valid array with P elements, /// where P = s.size()/src.size(). /// @param[out] produced must also point to a valid array with P elements. void computeInjectedProduced(const IncompPropertiesInterface& props, const std::vector<double>& s, const std::vector<double>& src, const double dt, double* injected, double* produced) { const int num_cells = src.size(); const int np = s.size()/src.size(); if (int(s.size()) != num_cells*np) { OPM_THROW(std::runtime_error, "Sizes of s and src vectors do not match."); } std::fill(injected, injected + np, 0.0); std::fill(produced, produced + np, 0.0); const double* visc = props.viscosity(); std::vector<double> mob(np); for (int c = 0; c < num_cells; ++c) { if (src[c] > 0.0) { injected[0] += src[c]*dt; } else if (src[c] < 0.0) { const double flux = -src[c]*dt; const double* sat = &s[np*c]; props.relperm(1, sat, &c, &mob[0], 0); double totmob = 0.0; for (int p = 0; p < np; ++p) { mob[p] /= visc[p]; totmob += mob[p]; } for (int p = 0; p < np; ++p) { produced[p] += (mob[p]/totmob)*flux; } } } }
void WellReport::push(const IncompPropertiesInterface& props, const Wells& wells, const std::vector<double>& saturation, const double time, const std::vector<double>& well_bhp, const std::vector<double>& well_perfrates) { int nw = well_bhp.size(); assert(nw == wells.number_of_wells); int np = props.numPhases(); const int max_np = 3; if (np > max_np) { OPM_THROW(std::runtime_error, "WellReport for now assumes #phases <= " << max_np); } const double* visc = props.viscosity(); std::vector<double> data_now; data_now.reserve(1 + 3*nw); data_now.push_back(time/unit::day); for (int w = 0; w < nw; ++w) { data_now.push_back(well_bhp[w]/(unit::barsa)); double well_rate_total = 0.0; double well_rate_water = 0.0; for (int perf = wells.well_connpos[w]; perf < wells.well_connpos[w + 1]; ++perf) { const double perf_rate = unit::convert::to(well_perfrates[perf], unit::cubic(unit::meter)/unit::day); well_rate_total += perf_rate; if (perf_rate > 0.0) { // Injection. well_rate_water += perf_rate*wells.comp_frac[0]; } else { // Production. const int cell = wells.well_cells[perf]; double mob[max_np]; props.relperm(1, &saturation[2*cell], &cell, mob, 0); double tmob = 0; for(int i = 0; i < np; ++i) { mob[i] /= visc[i]; tmob += mob[i]; } const double fracflow = mob[0]/tmob; well_rate_water += perf_rate*fracflow; } } data_now.push_back(well_rate_total); if (well_rate_total == 0.0) { data_now.push_back(0.0); } else { data_now.push_back(well_rate_water/well_rate_total); } } data_.push_back(data_now); }
/// @brief Computes injected and produced volumes of all phases, /// and injected and produced polymer mass. /// Note 1: assumes that only the first phase is injected. /// Note 2: assumes that transport has been done with an /// implicit method, i.e. that the current state /// gives the mobilities used for the preceding timestep. /// @param[in] props fluid and rock properties. /// @param[in] polyprops polymer properties /// @param[in] state state variables (pressure, fluxes etc.) /// @param[in] src if < 0: total reservoir volume outflow, /// if > 0: first phase reservoir volume inflow. /// @param[in] inj_c injected concentration by cell /// @param[in] dt timestep used /// @param[out] injected must point to a valid array with P elements, /// where P = s.size()/src.size(). /// @param[out] produced must also point to a valid array with P elements. /// @param[out] polyinj injected mass of polymer /// @param[out] polyprod produced mass of polymer void computeInjectedProduced(const IncompPropertiesInterface& props, const Opm::PolymerProperties& polyprops, const PolymerState& state, const std::vector<double>& transport_src, const std::vector<double>& inj_c, const double dt, double* injected, double* produced, double& polyinj, double& polyprod) { const int num_cells = transport_src.size(); if (props.numCells() != num_cells) { OPM_THROW(std::runtime_error, "Size of transport_src vector does not match number of cells in props."); } const int np = props.numPhases(); if (int(state.saturation().size()) != num_cells*np) { OPM_THROW(std::runtime_error, "Sizes of state vectors do not match number of cells."); } const std::vector<double>& s = state.saturation(); const std::vector<double>& c = state.getCellData( state.CONCENTRATION ); const std::vector<double>& cmax = state.getCellData( state.CMAX ); std::fill(injected, injected + np, 0.0); std::fill(produced, produced + np, 0.0); polyinj = 0.0; polyprod = 0.0; const double* visc = props.viscosity(); std::vector<double> kr_cell(np); double mob[2]; double mc; for (int cell = 0; cell < num_cells; ++cell) { if (transport_src[cell] > 0.0) { injected[0] += transport_src[cell]*dt; polyinj += transport_src[cell]*dt*inj_c[cell]; } else if (transport_src[cell] < 0.0) { const double flux = -transport_src[cell]*dt; const double* sat = &s[np*cell]; props.relperm(1, sat, &cell, &kr_cell[0], 0); polyprops.effectiveMobilities(c[cell], cmax[cell], visc, &kr_cell[0], mob); double totmob = mob[0] + mob[1]; for (int p = 0; p < np; ++p) { produced[p] += (mob[p]/totmob)*flux; } polyprops.computeMc(c[cell], mc); polyprod += (mob[0]/totmob)*flux*mc; } } }
/// @brief Computes total absorbed polymer mass over all grid cells. /// @param[in] props fluid and rock properties. /// @param[in] polyprops polymer properties /// @param[in] pv the pore volume by cell. /// @param[in] cmax max polymer concentration for cell /// @return total absorbed polymer mass. double computePolymerAdsorbed(const IncompPropertiesInterface& props, const Opm::PolymerProperties& polyprops, const std::vector<double>& pv, const std::vector<double>& cmax) { const int num_cells = pv.size(); const double rhor = polyprops.rockDensity(); const double* poro = props.porosity(); double abs_mass = 0.0; for (int cell = 0; cell < num_cells; ++cell) { double c_ads; polyprops.simpleAdsorption(cmax[cell], c_ads); abs_mass += c_ads*pv[cell]*((1.0 - poro[cell])/poro[cell])*rhor; } return abs_mass; }
void initStateBasic(const UnstructuredGrid& grid, const IncompPropertiesInterface& props, const parameter::ParameterGroup& param, const double gravity, State& state) { const int num_phases = props.numPhases(); if (num_phases != 2) { THROW("initStateTwophaseBasic(): currently handling only two-phase scenarios."); } state.init(grid, num_phases); const int num_cells = props.numCells(); // By default: initialise water saturation to minimum everywhere. std::vector<int> all_cells(num_cells); for (int i = 0; i < num_cells; ++i) { all_cells[i] = i; } state.setFirstSat(all_cells, props, State::MinSat); const bool convection_testcase = param.getDefault("convection_testcase", false); const bool segregation_testcase = param.getDefault("segregation_testcase", false); if (convection_testcase) { // Initialise water saturation to max in the 'left' part. std::vector<int> left_cells; left_cells.reserve(num_cells/2); const int *glob_cell = grid.global_cell; const int* cd = grid.cartdims; for (int cell = 0; cell < num_cells; ++cell) { const int gc = glob_cell == 0 ? cell : glob_cell[cell]; bool left = (gc % cd[0]) < cd[0]/2; if (left) { left_cells.push_back(cell); } } state.setFirstSat(left_cells, props, State::MaxSat); const double init_p = param.getDefault("ref_pressure", 100.0)*unit::barsa; std::fill(state.pressure().begin(), state.pressure().end(), init_p); } else if (segregation_testcase) { // Warn against error-prone usage. if (gravity == 0.0) { std::cout << "**** Warning: running gravity segregation scenario, but gravity is zero." << std::endl; } if (grid.cartdims[2] <= 1) { std::cout << "**** Warning: running gravity segregation scenario, which expects nz > 1." << std::endl; } // Initialise water saturation to max *above* water-oil contact. const double woc = param.get<double>("water_oil_contact"); initWaterOilContact(grid, props, woc, WaterAbove, state); // Initialise pressure to hydrostatic state. const double ref_p = param.getDefault("ref_pressure", 100.0)*unit::barsa; double dens[2] = { props.density()[1], props.density()[0] }; initHydrostaticPressure(grid, dens, woc, gravity, woc, ref_p, state); } else if (param.has("water_oil_contact")) { // Warn against error-prone usage. if (gravity == 0.0) { std::cout << "**** Warning: running gravity convection scenario, but gravity is zero." << std::endl; } if (grid.cartdims[2] <= 1) { std::cout << "**** Warning: running gravity convection scenario, which expects nz > 1." << std::endl; } // Initialise water saturation to max below water-oil contact. const double woc = param.get<double>("water_oil_contact"); initWaterOilContact(grid, props, woc, WaterBelow, state); // Initialise pressure to hydrostatic state. const double ref_p = param.getDefault("ref_pressure", 100.0)*unit::barsa; initHydrostaticPressure(grid, props.density(), woc, gravity, woc, ref_p, state); } else if (param.has("init_saturation")) { // Initialise water saturation to init_saturation parameter. const double init_saturation = param.get<double>("init_saturation"); for (int cell = 0; cell < num_cells; ++cell) { state.saturation()[2*cell] = init_saturation; state.saturation()[2*cell + 1] = 1.0 - init_saturation; } // Initialise pressure to hydrostatic state. const double ref_p = param.getDefault("ref_pressure", 100.0)*unit::barsa; const double rho = props.density()[0]*init_saturation + props.density()[1]*(1.0 - init_saturation); const double dens[2] = { rho, rho }; const double ref_z = grid.cell_centroids[0 + grid.dimensions - 1]; initHydrostaticPressure(grid, dens, ref_z, gravity, ref_z, ref_p, state); } else { // Use default: water saturation is minimum everywhere. // Initialise pressure to hydrostatic state. const double ref_p = param.getDefault("ref_pressure", 100.0)*unit::barsa; const double rho = props.density()[1]; const double dens[2] = { rho, rho }; const double ref_z = grid.cell_centroids[0 + grid.dimensions - 1]; initHydrostaticPressure(grid, dens, ref_z, gravity, ref_z, ref_p, state); } // Finally, init face pressures. initFacePressure(grid, state); }
SimulatorIncompTwophase::Impl::Impl(const parameter::ParameterGroup& param, const UnstructuredGrid& grid, const IncompPropertiesInterface& props, const RockCompressibility* rock_comp_props, WellsManager& wells_manager, const std::vector<double>& src, const FlowBoundaryConditions* bcs, LinearSolverInterface& linsolver, const double* gravity) : use_reorder_(param.getDefault("use_reorder", true)), use_segregation_split_(param.getDefault("use_segregation_split", false)), grid_(grid), props_(props), rock_comp_props_(rock_comp_props), wells_manager_(wells_manager), wells_(wells_manager.c_wells()), src_(src), bcs_(bcs), psolver_(grid, props, rock_comp_props, linsolver, param.getDefault("nl_pressure_residual_tolerance", 0.0), param.getDefault("nl_pressure_change_tolerance", 1.0), param.getDefault("nl_pressure_maxiter", 10), gravity, wells_manager.c_wells(), src, bcs) { // Initialize transport solver. if (use_reorder_) { tsolver_.reset(new Opm::TransportSolverTwophaseReorder(grid, props, use_segregation_split_ ? gravity : NULL, param.getDefault("nl_tolerance", 1e-9), param.getDefault("nl_maxiter", 30))); } else { if (rock_comp_props && rock_comp_props->isActive()) { OPM_THROW(std::runtime_error, "The implicit pressure solver cannot handle rock compressibility."); } if (use_segregation_split_) { OPM_THROW(std::runtime_error, "The implicit pressure solver is not set up to use segregation splitting."); } std::vector<double> porevol; computePorevolume(grid, props.porosity(), porevol); tsolver_.reset(new Opm::TransportSolverTwophaseImplicit(grid, props, porevol, gravity, psolver_.getHalfTrans(), param)); } // For output. log_ = param.getDefault("quiet", false) ? &Opm::null_stream : &std::cout; output_ = param.getDefault("output", true); if (output_) { output_vtk_ = param.getDefault("output_vtk", true); output_dir_ = param.getDefault("output_dir", std::string("output")); // Ensure that output dir exists boost::filesystem::path fpath(output_dir_); try { create_directories(fpath); } catch (...) { OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath); } output_interval_ = param.getDefault("output_interval", 1); } // Well control related init. check_well_controls_ = param.getDefault("check_well_controls", false); max_well_control_iterations_ = param.getDefault("max_well_control_iterations", 10); // Transport related init. num_transport_substeps_ = param.getDefault("num_transport_substeps", 1); // Misc init. const int num_cells = grid.number_of_cells; allcells_.resize(num_cells); for (int cell = 0; cell < num_cells; ++cell) { allcells_[cell] = cell; } }