/// Gas viscosity. /// \param[in] pg Array of n gas pressure values. /// \param[in] rv Array of n vapor oil/gas ratio /// \param[in] cond Array of n objects, each specifying which phases are present with non-zero saturation in a cell. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n viscosity values. ADB BlackoilPropsAd::muGas(const ADB& pg, const ADB& rv, const std::vector<PhasePresence>& cond, const Cells& cells) const { #if 1 return ADB::constant(muGas(pg.value(), rv.value(),cond,cells), pg.blockPattern()); #else if (!pu_.phase_used[Gas]) { OPM_THROW(std::runtime_error, "Cannot call muGas(): gas phase not present."); } const int n = cells.size(); assert(pg.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); if (pu_.phase_used[Oil]) { // Faking a z with the right ratio: // rv = zo/zg z.col(pu_.phase_pos[Oil]) = rv; z.col(pu_.phase_pos[Gas]) = V::Ones(n, 1); } Block mu(n, np); Block dmu(n, np); props_.viscosity(n, pg.value().data(), z.data(), cells.data(), mu.data(), dmu.data()); ADB::M dmu_diag = spdiag(dmu.col(pu_.phase_pos[Gas])); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmu_diag * pg.derivative()[block]; } return ADB::function(mu.col(pu_.phase_pos[Gas]), jacs); #endif }
/// Water viscosity. /// \param[in] pw Array of n water pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n viscosity values. ADB BlackoilPropsAd::muWat(const ADB& pw, const Cells& cells) const { #if 1 return ADB::constant(muWat(pw.value(), cells), pw.blockPattern()); #else if (!pu_.phase_used[Water]) { OPM_THROW(std::runtime_error, "Cannot call muWat(): water phase not present."); } const int n = cells.size(); assert(pw.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); Block mu(n, np); Block dmu(n, np); props_.viscosity(n, pw.value().data(), z.data(), cells.data(), mu.data(), dmu.data()); ADB::M dmu_diag = spdiag(dmu.col(pu_.phase_pos[Water])); const int num_blocks = pw.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmu_diag * pw.derivative()[block]; } return ADB::function(mu.col(pu_.phase_pos[Water]), jacs); #endif }
ADB SolventPropsAdFromDeck::muSolvent(const ADB& pg, const Cells& cells) const { const int n = cells.size(); assert(pg.value().size() == n); V mu(n); V dmudp(n); for (int i = 0; i < n; ++i) { const double& pg_i = pg.value()[i]; int regionIdx = cellPvtRegionIdx_[cells[i]]; double tempInvB = b_[regionIdx](pg_i); double tempInvBmu = inverseBmu_[regionIdx](pg_i); mu[i] = tempInvB / tempInvBmu; dmudp[i] = (tempInvBmu * b_[regionIdx].derivative(pg_i) - tempInvB * inverseBmu_[regionIdx].derivative(pg_i)) / (tempInvBmu * tempInvBmu); } ADB::M dmudp_diag(dmudp.matrix().asDiagonal()); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmudp_diag * pg.derivative()[block]; } return ADB::function(std::move(mu), std::move(jacs)); }
/// Oil formation volume factor. /// \param[in] po Array of n oil pressure values. /// \param[in] rs Array of n gas solution factor values. /// \param[in] cond Array of n taxonomies classifying fluid condition. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAd::bOil(const ADB& po, const ADB& rs, const std::vector<PhasePresence>& /*cond*/, const Cells& cells) const { if (!pu_.phase_used[Oil]) { OPM_THROW(std::runtime_error, "Cannot call muOil(): oil phase not present."); } const int n = cells.size(); assert(po.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); if (pu_.phase_used[Gas]) { // Faking a z with the right ratio: // rs = zg/zo z.col(pu_.phase_pos[Oil]) = V::Ones(n, 1); z.col(pu_.phase_pos[Gas]) = rs.value(); } Block matrix(n, np*np); Block dmatrix(n, np*np); props_.matrix(n, po.value().data(), z.data(), cells.data(), matrix.data(), dmatrix.data()); const int phase_ind = pu_.phase_pos[Oil]; const int column = phase_ind*np + phase_ind; // Index of our sought diagonal column. ADB::M db_diag = spdiag(dmatrix.col(column)); const int num_blocks = po.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { // For now, we deliberately ignore the derivative with respect to rs, // since the BlackoilPropertiesInterface class does not evaluate it. // We would add to the next line: + db_drs_diag * rs.derivative()[block] jacs[block] = db_diag * po.derivative()[block]; } return ADB::function(matrix.col(column), jacs); }
/// Gas viscosity. /// \param[in] pg Array of n gas pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n viscosity values. ADB BlackoilPropsAd::muGas(const ADB& pg, const Cells& cells) const { #if 1 return ADB::constant(muGas(pg.value(), cells), pg.blockPattern()); #else if (!pu_.phase_used[Gas]) { THROW("Cannot call muGas(): gas phase not present."); } const int n = cells.size(); ASSERT(pg.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); Block mu(n, np); Block dmu(n, np); props_.viscosity(n, pg.value().data(), z.data(), cells.data(), mu.data(), dmu.data()); ADB::M dmu_diag = spdiag(dmu.col(pu_.phase_pos[Gas])); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmu_diag * pg.derivative()[block]; } return ADB::function(mu.col(pu_.phase_pos[Gas]), jacs); #endif }
/// Gas formation volume factor. /// \param[in] pg Array of n gas pressure values. /// \param[in] rv Array of n vapor oil/gas ratio /// \param[in] cond Array of n objects, each specifying which phases are present with non-zero saturation in a cell. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAd::bGas(const ADB& pg, const ADB& rv, const std::vector<PhasePresence>& /*cond*/, const Cells& cells) const { if (!pu_.phase_used[Gas]) { OPM_THROW(std::runtime_error, "Cannot call muGas(): gas phase not present."); } const int n = cells.size(); assert(pg.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); if (pu_.phase_used[Oil]) { // Faking a z with the right ratio: // rv = zo/zg z.col(pu_.phase_pos[Oil]) = rv.value(); z.col(pu_.phase_pos[Gas]) = V::Ones(n, 1); } Block matrix(n, np*np); Block dmatrix(n, np*np); props_.matrix(n, pg.value().data(), z.data(), cells.data(), matrix.data(), dmatrix.data()); const int phase_ind = pu_.phase_pos[Gas]; const int column = phase_ind*np + phase_ind; // Index of our sought diagonal column. ADB::M db_diag = spdiag(dmatrix.col(column)); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = db_diag * pg.derivative()[block]; } return ADB::function(matrix.col(column), jacs); }
/// Gas formation volume factor. /// \param[in] pg Array of n gas pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAdFromDeck::bGas(const ADB& pg, const Cells& cells) const { if (!phase_usage_.phase_used[Gas]) { OPM_THROW(std::runtime_error, "Cannot call muGas(): gas phase not present."); } const int n = cells.size(); assert(pg.size() == n); V b(n); V dbdp(n); V dbdr(n); const double* rs = 0; props_[phase_usage_.phase_pos[Gas]]->b(n, pg.value().data(), rs, b.data(), dbdp.data(), dbdr.data()); ADB::M dbdp_diag = spdiag(dbdp); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dbdp_diag * pg.derivative()[block]; } return ADB::function(b, jacs); }
/// Oil viscosity. /// \param[in] po Array of n oil pressure values. /// \param[in] rs Array of n gas solution factor values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n viscosity values. ADB BlackoilPropsAd::muOil(const ADB& po, const ADB& rs, const Cells& cells) const { #if 1 return ADB::constant(muOil(po.value(), rs.value(), cells), po.blockPattern()); #else if (!pu_.phase_used[Oil]) { THROW("Cannot call muOil(): oil phase not present."); } const int n = cells.size(); ASSERT(po.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); if (pu_.phase_used[Gas]) { // Faking a z with the right ratio: // rs = zg/zo z.col(pu_.phase_pos[Oil]) = V::Ones(n, 1); z.col(pu_.phase_pos[Gas]) = rs.value(); } Block mu(n, np); Block dmu(n, np); props_.viscosity(n, po.value().data(), z.data(), cells.data(), mu.data(), dmu.data()); ADB::M dmu_diag = spdiag(dmu.col(pu_.phase_pos[Oil])); const int num_blocks = po.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { // For now, we deliberately ignore the derivative with respect to rs, // since the BlackoilPropertiesInterface class does not evaluate it. // We would add to the next line: + dmu_drs_diag * rs.derivative()[block] jacs[block] = dmu_diag * po.derivative()[block]; } return ADB::function(mu.col(pu_.phase_pos[Oil]), jacs); #endif }
/// Oil formation volume factor. /// \param[in] po Array of n oil pressure values. /// \param[in] rs Array of n gas solution factor values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAdFromDeck::bOil(const ADB& po, const ADB& rs, const Cells& cells) const { if (!phase_usage_.phase_used[Oil]) { OPM_THROW(std::runtime_error, "Cannot call muOil(): oil phase not present."); } const int n = cells.size(); assert(po.size() == n); V b(n); V dbdp(n); V dbdr(n); props_[phase_usage_.phase_pos[Oil]]->b(n, po.value().data(), rs.value().data(), b.data(), dbdp.data(), dbdr.data()); ADB::M dbdp_diag = spdiag(dbdp); ADB::M dbdr_diag = spdiag(dbdr); const int num_blocks = po.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dbdp_diag * po.derivative()[block] + dbdr_diag * rs.derivative()[block]; } return ADB::function(b, jacs); }
ADB PolymerPropsAd::polymerWaterVelocityRatio(const ADB& c) const { const int nc = c.size(); V mc(nc); V dmc(nc); for (int i = 0; i < nc; ++i) { double m = 0; double dm = 0; polymer_props_.computeMcWithDer(c.value()(i), m, dm); mc(i) = m; dmc(i) = dm; } ADB::M dmc_diag(dmc.matrix().asDiagonal()); const int num_blocks = c.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmc_diag * c.derivative()[block]; } return ADB::function(std::move(mc), std::move(jacs)); }
std::vector<ADB> BlackoilPropsAd::capPress(const ADB& sw, const ADB& so, const ADB& sg, const Cells& cells) const { const int numCells = cells.size(); const int numActivePhases = numPhases(); const int numBlocks = so.numBlocks(); Block activeSat(numCells, numActivePhases); if (pu_.phase_used[Water]) { assert(sw.value().size() == numCells); activeSat.col(pu_.phase_pos[Water]) = sw.value(); } if (pu_.phase_used[Oil]) { assert(so.value().size() == numCells); activeSat.col(pu_.phase_pos[Oil]) = so.value(); } else { OPM_THROW(std::runtime_error, "BlackoilPropsAdFromDeck::relperm() assumes oil phase is active."); } if (pu_.phase_used[Gas]) { assert(sg.value().size() == numCells); activeSat.col(pu_.phase_pos[Gas]) = sg.value(); } Block pc(numCells, numActivePhases); Block dpc(numCells, numActivePhases*numActivePhases); props_.capPress(numCells, activeSat.data(), cells.data(), pc.data(), dpc.data()); std::vector<ADB> adbCapPressures; adbCapPressures.reserve(3); const ADB* s[3] = { &sw, &so, &sg }; for (int phase1 = 0; phase1 < 3; ++phase1) { if (pu_.phase_used[phase1]) { const int phase1_pos = pu_.phase_pos[phase1]; std::vector<ADB::M> jacs(numBlocks); for (int block = 0; block < numBlocks; ++block) { jacs[block] = ADB::M(numCells, s[phase1]->derivative()[block].cols()); } for (int phase2 = 0; phase2 < 3; ++phase2) { if (!pu_.phase_used[phase2]) continue; const int phase2_pos = pu_.phase_pos[phase2]; // Assemble dpc1/ds2. const int column = phase1_pos + numActivePhases*phase2_pos; // Recall: Fortran ordering from props_.relperm() ADB::M dpc1_ds2_diag = spdiag(dpc.col(column)); for (int block = 0; block < numBlocks; ++block) { jacs[block] += dpc1_ds2_diag * s[phase2]->derivative()[block]; } } adbCapPressures.emplace_back(ADB::function(pc.col(phase1_pos), jacs)); } else { adbCapPressures.emplace_back(ADB::null()); } } return adbCapPressures; }
/// Relative permeabilities for all phases. /// \param[in] sw Array of n water saturation values. /// \param[in] so Array of n oil saturation values. /// \param[in] sg Array of n gas saturation values. /// \param[in] cells Array of n cell indices to be associated with the saturation values. /// \return An std::vector with 3 elements, each an array of n relperm values, /// containing krw, kro, krg. Use PhaseIndex for indexing into the result. std::vector<ADB> BlackoilPropsAd::relperm(const ADB& sw, const ADB& so, const ADB& sg, const Cells& cells) const { const int n = cells.size(); const int np = props_.numPhases(); Block s_all(n, np); if (pu_.phase_used[Water]) { assert(sw.value().size() == n); s_all.col(pu_.phase_pos[Water]) = sw.value(); } if (pu_.phase_used[Oil]) { assert(so.value().size() == n); s_all.col(pu_.phase_pos[Oil]) = so.value(); } else { OPM_THROW(std::runtime_error, "BlackoilPropsAd::relperm() assumes oil phase is active."); } if (pu_.phase_used[Gas]) { assert(sg.value().size() == n); s_all.col(pu_.phase_pos[Gas]) = sg.value(); } Block kr(n, np); Block dkr(n, np*np); props_.relperm(n, s_all.data(), cells.data(), kr.data(), dkr.data()); const int num_blocks = so.numBlocks(); std::vector<ADB> relperms; relperms.reserve(3); typedef const ADB* ADBPtr; ADBPtr s[3] = { &sw, &so, &sg }; for (int phase1 = 0; phase1 < 3; ++phase1) { if (pu_.phase_used[phase1]) { const int phase1_pos = pu_.phase_pos[phase1]; std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = ADB::M(n, s[phase1]->derivative()[block].cols()); } for (int phase2 = 0; phase2 < 3; ++phase2) { if (!pu_.phase_used[phase2]) { continue; } const int phase2_pos = pu_.phase_pos[phase2]; // Assemble dkr1/ds2. const int column = phase1_pos + np*phase2_pos; // Recall: Fortran ordering from props_.relperm() ADB::M dkr1_ds2_diag = spdiag(dkr.col(column)); for (int block = 0; block < num_blocks; ++block) { jacs[block] += dkr1_ds2_diag * s[phase2]->derivative()[block]; } } relperms.emplace_back(ADB::function(kr.col(phase1_pos), jacs)); } else { relperms.emplace_back(ADB::null()); } } return relperms; }
/// Bubble point curve for Rs as function of oil pressure. /// \param[in] po Array of n oil pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n bubble point values for Rs. ADB BlackoilPropsAdFromDeck::rsMax(const ADB& po, const Cells& cells) const { if (!phase_usage_.phase_used[Oil]) { OPM_THROW(std::runtime_error, "Cannot call rsMax(): oil phase not present."); } const int n = cells.size(); assert(po.size() == n); V rbub(n); V drbubdp(n); props_[Oil]->rbub(n, po.value().data(), rbub.data(), drbubdp.data()); ADB::M drbubdp_diag = spdiag(drbubdp); const int num_blocks = po.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = drbubdp_diag * po.derivative()[block]; } return ADB::function(rbub, jacs); }
ADB PolymerPropsAd::effectiveInvWaterVisc(const ADB& c, const double* visc) const { const int nc = c.size(); V inv_mu_w_eff(nc); V dinv_mu_w_eff(nc); for (int i = 0; i < nc; ++i) { double im = 0, dim = 0; polymer_props_.effectiveInvViscWithDer(c.value()(i), visc, im, dim); inv_mu_w_eff(i) = im; dinv_mu_w_eff(i) = dim; } ADB::M dim_diag = spdiag(dinv_mu_w_eff); const int num_blocks = c.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dim_diag * c.derivative()[block]; } return ADB::function(std::move(inv_mu_w_eff), std::move(jacs)); }
ADB PolymerPropsAd::viscMult(const ADB& c) const { const int nc = c.size(); V visc_mult(nc); V dvisc_mult(nc); for (int i = 0; i < nc; ++i) { double im = 0, dim = 0; im = polymer_props_.viscMultWithDer(c.value()(i), &dim); visc_mult(i) = im; dvisc_mult(i) = dim; } ADB::M dim_diag(dvisc_mult.matrix().asDiagonal()); const int num_blocks = c.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dim_diag * c.derivative()[block]; } return ADB::function(std::move(visc_mult), std::move(jacs)); }
ADB SolventPropsAdFromDeck::makeADBfromTables(const ADB& X_AD, const Cells& cells, const std::vector<int>& regionIdx, const std::vector<NonuniformTableLinear<double>>& tables) const { const int n = cells.size(); assert(X_AD.value().size() == n); V x(n); V dx(n); for (int i = 0; i < n; ++i) { const double& X_i = X_AD.value()[i]; x[i] = tables[regionIdx[cells[i]]](X_i); dx[i] = tables[regionIdx[cells[i]]].derivative(X_i); } ADB::M dx_diag(dx.matrix().asDiagonal()); const int num_blocks = X_AD.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { fastSparseProduct(dx_diag, X_AD.derivative()[block], jacs[block]); } return ADB::function(std::move(x), std::move(jacs)); }
/// Gas formation volume factor. /// \param[in] pg Array of n gas pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAd::bGas(const ADB& pg, const Cells& cells) const { if (!pu_.phase_used[Gas]) { THROW("Cannot call muGas(): gas phase not present."); } const int n = cells.size(); ASSERT(pg.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); Block matrix(n, np*np); Block dmatrix(n, np*np); props_.matrix(n, pg.value().data(), z.data(), cells.data(), matrix.data(), dmatrix.data()); const int phase_ind = pu_.phase_pos[Gas]; const int column = phase_ind*np + phase_ind; // Index of our sought diagonal column. ADB::M db_diag = spdiag(dmatrix.col(column)); const int num_blocks = pg.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = db_diag * pg.derivative()[block]; } return ADB::function(matrix.col(column), jacs); }
/// Water formation volume factor. /// \param[in] pw Array of n water pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n formation volume factor values. ADB BlackoilPropsAd::bWat(const ADB& pw, const Cells& cells) const { if (!pu_.phase_used[Water]) { OPM_THROW(std::runtime_error, "Cannot call muWat(): water phase not present."); } const int n = cells.size(); assert(pw.value().size() == n); const int np = props_.numPhases(); Block z = Block::Zero(n, np); Block matrix(n, np*np); Block dmatrix(n, np*np); props_.matrix(n, pw.value().data(), z.data(), cells.data(), matrix.data(), dmatrix.data()); const int phase_ind = pu_.phase_pos[Water]; const int column = phase_ind*np + phase_ind; // Index of our sought diagonal column. ADB::M db_diag = spdiag(dmatrix.col(column)); const int num_blocks = pw.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = db_diag * pw.derivative()[block]; } return ADB::function(matrix.col(column), jacs); }
/// Water viscosity. /// \param[in] pw Array of n water pressure values. /// \param[in] cells Array of n cell indices to be associated with the pressure values. /// \return Array of n viscosity values. ADB BlackoilPropsAdFromDeck::muWat(const ADB& pw, const Cells& cells) const { if (!phase_usage_.phase_used[Water]) { OPM_THROW(std::runtime_error, "Cannot call muWat(): water phase not present."); } const int n = cells.size(); assert(pw.size() == n); V mu(n); V dmudp(n); V dmudr(n); const double* rs = 0; props_[phase_usage_.phase_pos[Water]]->mu(n, pw.value().data(), rs, mu.data(), dmudp.data(), dmudr.data()); ADB::M dmudp_diag = spdiag(dmudp); const int num_blocks = pw.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dmudp_diag * pw.derivative()[block]; } return ADB::function(mu, jacs); }
ADB PolymerPropsAd::effectiveInvPolymerVisc(const ADB& c, const V& mu_w) const { assert(c.size() == mu_w.size()); const int nc = c.size(); V inv_mu_p_eff(nc); V dinv_mu_p_eff(nc); for (int i = 0; i < nc; ++i) { double im = 0; double dim = 0; polymer_props_.effectiveInvPolyViscWithDer(c.value()(i), mu_w(i), im, dim); inv_mu_p_eff(i) = im; dinv_mu_p_eff(i) = dim; } ADB::M dim_diag(dinv_mu_p_eff.matrix().asDiagonal()); const int num_blocks = c.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dim_diag * c.derivative()[block]; } return ADB::function(std::move(inv_mu_p_eff), std::move(jacs)); }
ADB PolymerPropsAd::adsorption(const ADB& c, const ADB& cmax_cells) const { const int nc = c.value().size(); V ads(nc); V dads(nc); for (int i = 0; i < nc; ++i) { double c_ads = 0; double dc_ads = 0; polymer_props_.adsorptionWithDer(c.value()(i), cmax_cells.value()(i), c_ads, dc_ads); ads(i) = c_ads; dads(i) = dc_ads; } ADB::M dads_diag(dads.matrix().asDiagonal()); int num_blocks = c.numBlocks(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { jacs[block] = dads_diag * c.derivative()[block]; } return ADB::function(std::move(ads), std::move(jacs)); }
VFPProdProperties::ADB VFPProdProperties::bhp(const std::vector<int>& table_id, const ADB& aqua, const ADB& liquid, const ADB& vapour, const ADB& thp_arg, const ADB& alq) const { const int nw = thp_arg.size(); std::vector<int> block_pattern = detail::commonBlockPattern(aqua, liquid, vapour, thp_arg, alq); assert(static_cast<int>(table_id.size()) == nw); assert(aqua.size() == nw); assert(liquid.size() == nw); assert(vapour.size() == nw); assert(thp_arg.size() == nw); assert(alq.size() == nw); //Allocate data for bhp's and partial derivatives ADB::V value = ADB::V::Zero(nw); ADB::V dthp = ADB::V::Zero(nw); ADB::V dwfr = ADB::V::Zero(nw); ADB::V dgfr = ADB::V::Zero(nw); ADB::V dalq = ADB::V::Zero(nw); ADB::V dflo = ADB::V::Zero(nw); //Get the table for each well std::vector<const VFPProdTable*> well_tables(nw, nullptr); for (int i=0; i<nw; ++i) { if (table_id[i] >= 0) { well_tables[i] = detail::getTable(m_tables, table_id[i]); } } //Get the right FLO/GFR/WFR variable for each well as a single ADB const ADB flo = detail::combineADBVars<VFPProdTable::FLO_TYPE>(well_tables, aqua, liquid, vapour); const ADB wfr = detail::combineADBVars<VFPProdTable::WFR_TYPE>(well_tables, aqua, liquid, vapour); const ADB gfr = detail::combineADBVars<VFPProdTable::GFR_TYPE>(well_tables, aqua, liquid, vapour); //Compute the BHP for each well independently for (int i=0; i<nw; ++i) { const VFPProdTable* table = well_tables[i]; if (table != nullptr) { //First, find the values to interpolate between //Value of FLO is negative in OPM for producers, but positive in VFP table auto flo_i = detail::findInterpData(-flo.value()[i], table->getFloAxis()); auto thp_i = detail::findInterpData( thp_arg.value()[i], table->getTHPAxis()); auto wfr_i = detail::findInterpData( wfr.value()[i], table->getWFRAxis()); auto gfr_i = detail::findInterpData( gfr.value()[i], table->getGFRAxis()); auto alq_i = detail::findInterpData( alq.value()[i], table->getALQAxis()); detail::VFPEvaluation bhp_val = detail::interpolate(table->getTable(), flo_i, thp_i, wfr_i, gfr_i, alq_i); value[i] = bhp_val.value; dthp[i] = bhp_val.dthp; dwfr[i] = bhp_val.dwfr; dgfr[i] = bhp_val.dgfr; dalq[i] = bhp_val.dalq; dflo[i] = bhp_val.dflo; } else { value[i] = -1e100; //Signal that this value has not been calculated properly, due to "missing" table } } //Create diagonal matrices from ADB::Vs ADB::M dthp_diag(dthp.matrix().asDiagonal()); ADB::M dwfr_diag(dwfr.matrix().asDiagonal()); ADB::M dgfr_diag(dgfr.matrix().asDiagonal()); ADB::M dalq_diag(dalq.matrix().asDiagonal()); ADB::M dflo_diag(dflo.matrix().asDiagonal()); //Calculate the Jacobians const int num_blocks = block_pattern.size(); std::vector<ADB::M> jacs(num_blocks); for (int block = 0; block < num_blocks; ++block) { //Could have used fastSparseProduct and temporary variables //but may not save too much on that. jacs[block] = ADB::M(nw, block_pattern[block]); if (!thp_arg.derivative().empty()) { jacs[block] += dthp_diag * thp_arg.derivative()[block]; } if (!wfr.derivative().empty()) { jacs[block] += dwfr_diag * wfr.derivative()[block]; } if (!gfr.derivative().empty()) { jacs[block] += dgfr_diag * gfr.derivative()[block]; } if (!alq.derivative().empty()) { jacs[block] += dalq_diag * alq.derivative()[block]; } if (!flo.derivative().empty()) { jacs[block] -= dflo_diag * flo.derivative()[block]; } } ADB retval = ADB::function(std::move(value), std::move(jacs)); return retval; }