/// Constructor SinglePvtDeadSpline::SinglePvtDeadSpline(const table_t& pvd_table, const int samples) { const int region_number = 0; if (pvd_table.size() != 1) { THROW("More than one PVT-region"); } // Copy data const int sz = pvd_table[region_number][0].size(); std::vector<double> press(sz); std::vector<double> B_inv(sz); std::vector<double> visc(sz); for (int i = 0; i < sz; ++i) { press[i] = pvd_table[region_number][0][i]; B_inv[i] = 1.0 / pvd_table[region_number][1][i]; visc[i] = pvd_table[region_number][2][i]; } buildUniformMonotoneTable(press, B_inv, samples, one_over_B_); buildUniformMonotoneTable(press, visc, samples, viscosity_); // Dumping the created tables. // static int count = 0; // std::ofstream os((std::string("dump-") + boost::lexical_cast<std::string>(count++)).c_str()); // os.precision(15); // os << "1/B\n\n" << one_over_B_ // << "\n\nvisc\n\n" << viscosity_ << std::endl; }
//------------------------------------------------------------------------ // Member functions //------------------------------------------------------------------------- /// Constructor SinglePvtDead::SinglePvtDead(const table_t& pvd_table) { const int region_number = 0; if (pvd_table.size() != 1) { OPM_THROW(std::runtime_error, "More than one PVT-region"); } // Copy data const int sz = pvd_table[region_number][0].size(); std::vector<double> press(sz); std::vector<double> b(sz); std::vector<double> visc(sz); for (int i = 0; i < sz; ++i) { press[i] = pvd_table[region_number][0][i]; b[i] = 1.0 / pvd_table[region_number][1][i]; visc[i] = pvd_table[region_number][2][i]; } b_ = NonuniformTableLinear<double>(press, b); viscosity_ = NonuniformTableLinear<double>(press, visc); // Dumping the created tables. // static int count = 0; // std::ofstream os((std::string("dump-") + boost::lexical_cast<std::string>(count++)).c_str()); // os.precision(15); // os << "1/B\n\n" << one_over_B_ // << "\n\nvisc\n\n" << viscosity_ << std::endl; }
/// @brief Computes injected and produced surface 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. /// Note 3: Gives surface volume values, not reservoir volumes /// (as the incompressible version of the function does). /// Also, assumes that transport_src is given in surface volumes /// for injector terms! /// @param[in] props fluid and rock properties. /// @param[in] state state variables (pressure, sat, surfvol) /// @param[in] transport_src if < 0: total resv outflow, if > 0: first phase surfv 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 BlackoilPropertiesInterface& props, const BlackoilState& state, const std::vector<double>& transport_src, const double dt, double* injected, double* produced) { 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>& press = state.pressure(); const std::vector<double>& temp = state.temperature(); const std::vector<double>& s = state.saturation(); const std::vector<double>& z = state.surfacevol(); std::fill(injected, injected + np, 0.0); std::fill(produced, produced + np, 0.0); std::vector<double> visc(np); std::vector<double> mob(np); std::vector<double> A(np*np); std::vector<double> prod_resv_phase(np); std::vector<double> prod_surfvol(np); for (int c = 0; c < num_cells; ++c) { if (transport_src[c] > 0.0) { // Inflowing transport source is a surface volume flux // for the first phase. injected[0] += transport_src[c]*dt; } else if (transport_src[c] < 0.0) { // Outflowing transport source is a total reservoir // volume flux. const double flux = -transport_src[c]*dt; const double* sat = &s[np*c]; props.relperm(1, sat, &c, &mob[0], 0); props.viscosity(1, &press[c], &temp[c], &z[np*c], &c, &visc[0], 0); props.matrix(1, &press[c], &temp[c], &z[np*c], &c, &A[0], 0); double totmob = 0.0; for (int p = 0; p < np; ++p) { mob[p] /= visc[p]; totmob += mob[p]; } std::fill(prod_surfvol.begin(), prod_surfvol.end(), 0.0); for (int p = 0; p < np; ++p) { prod_resv_phase[p] = (mob[p]/totmob)*flux; for (int q = 0; q < np; ++q) { prod_surfvol[q] += prod_resv_phase[p]*A[q + np*p]; } } for (int p = 0; p < np; ++p) { produced[p] += prod_surfvol[p]; } } } }
/// @brief Computes injected and produced volumes of all phases, /// and injected and produced polymer mass - in the compressible case. /// 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] transport_src if < 0: total reservoir volume outflow, /// if > 0: first phase *surface 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()/transport_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 BlackoilPropertiesInterface& props, const Opm::PolymerProperties& polyprops, const PolymerBlackoilState& 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>& press = state.pressure(); const std::vector<double>& temp = state.temperature(); const std::vector<double>& s = state.saturation(); const std::vector<double>& z = state.surfacevol(); 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; std::vector<double> visc(np); std::vector<double> kr_cell(np); std::vector<double> mob(np); std::vector<double> A(np*np); std::vector<double> prod_resv_phase(np); std::vector<double> prod_surfvol(np); double mc; for (int cell = 0; cell < num_cells; ++cell) { if (transport_src[cell] > 0.0) { // Inflowing transport source is a surface volume flux // for the first phase. injected[0] += transport_src[cell]*dt; polyinj += transport_src[cell]*dt*inj_c[cell]; } else if (transport_src[cell] < 0.0) { // Outflowing transport source is a total reservoir // volume flux. const double flux = -transport_src[cell]*dt; const double* sat = &s[np*cell]; props.relperm(1, sat, &cell, &kr_cell[0], 0); props.viscosity(1, &press[cell], &temp[cell], &z[np*cell], &cell, &visc[0], 0); props.matrix(1, &press[cell], &temp[cell], &z[np*cell], &cell, &A[0], 0); polyprops.effectiveMobilities(c[cell], cmax[cell], &visc[0], &kr_cell[0], &mob[0]); double totmob = 0.0; for (int p = 0; p < np; ++p) { totmob += mob[p]; } std::fill(prod_surfvol.begin(), prod_surfvol.end(), 0.0); for (int p = 0; p < np; ++p) { prod_resv_phase[p] = (mob[p]/totmob)*flux; for (int q = 0; q < np; ++q) { prod_surfvol[q] += prod_resv_phase[p]*A[q + np*p]; } } for (int p = 0; p < np; ++p) { produced[p] += prod_surfvol[p]; } polyprops.computeMc(c[cell], mc); polyprod += produced[0]*mc; } } }