예제 #1
0
    /// @brief Computes total mobility and omega for a set of s/c values.
    /// @param[in]  props     rock and fluid properties
    /// @param[in]  polyprops polymer properties
    /// @param[in]  cells     cells with which the saturation values are associated
    /// @param[in]  s         saturation values (for all phases)
    /// @param[in]  c         polymer concentration
    /// @param[out] totmob    total mobility
    /// @param[out] omega     mobility-weighted (or fractional-flow weighted)
    ///                       fluid densities.
    void computeTotalMobilityOmega(const Opm::IncompPropertiesInterface& props,
				   const Opm::PolymerProperties& polyprops,
				   const std::vector<int>& cells,
				   const std::vector<double>& s,
				   const std::vector<double>& c,
                                   const std::vector<double>& cmax,
				   std::vector<double>& totmob,
				   std::vector<double>& omega)
    {
	int num_cells = cells.size();
	int num_phases = props.numPhases();
	totmob.resize(num_cells);
	omega.resize(num_cells);
	assert(int(s.size()) == num_cells*num_phases);
	std::vector<double> kr(num_cells*num_phases);
	props.relperm(num_cells, &s[0], &cells[0], &kr[0], 0);
	const double* visc = props.viscosity();
	const double* rho = props.density();
        double mob[2]; // here we assume num_phases=2
	for (int cell = 0; cell < num_cells; ++cell) {
            double*  kr_cell = &kr[2*cell];
            polyprops.effectiveMobilities(c[cell], cmax[cell], visc, kr_cell,
                                          mob);
            totmob[cell] = mob[0] + mob[1];
            omega[cell] = rho[0]*mob[0]/totmob[cell] + rho[1]*mob[1]/totmob[cell];
        }
    }
예제 #2
0
    /// Computes the fractional flow for each cell in the cells argument
    /// @param[in]  props            rock and fluid properties
    /// @param[in]  polyprops        polymer properties
    /// @param[in]  cells            cells with which the saturation values are associated
    /// @param[in]  s                saturation values (for all phases)
    /// @param[in]  c                concentration values
    /// @param[in]  cmax             max polymer concentration experienced by cell
    /// @param[out] fractional_flow  the fractional flow for each phase for each cell.
    void computeFractionalFlow(const Opm::IncompPropertiesInterface& props,
                               const Opm::PolymerProperties& polyprops,
                               const std::vector<int>& cells,
                               const std::vector<double>& s,
                               const std::vector<double>& c,
                               const std::vector<double>& cmax,
                               std::vector<double>& fractional_flows)
    {
	int num_cells = cells.size();
	int num_phases = props.numPhases();
        if (num_phases != 2) {
            OPM_THROW(std::runtime_error, "computeFractionalFlow() assumes 2 phases.");
        }
	fractional_flows.resize(num_cells*num_phases);
	assert(int(s.size()) == num_cells*num_phases);
	std::vector<double> kr(num_cells*num_phases);
	props.relperm(num_cells, &s[0], &cells[0], &kr[0], 0);
	const double* visc = props.viscosity();
        double mob[2]; // here we assume num_phases=2
	for (int cell = 0; cell < num_cells; ++cell) {
            double* kr_cell = &kr[2*cell];
            polyprops.effectiveMobilities(c[cell], cmax[cell], visc, kr_cell, mob);
            fractional_flows[2*cell]     = mob[0] / (mob[0] + mob[1]);
            fractional_flows[2*cell + 1] = mob[1] / (mob[0] + mob[1]);
        }
    }
예제 #3
0
    /// @brief Computes total mobility for a set of s/c values.
    /// @param[in]  props     rock and fluid properties
    /// @param[in]  polyprops polymer properties
    /// @param[in]  cells     cells with which the saturation values are associated
    /// @param[in]  s         saturation values (for all phases)
    /// @param[in]  c         polymer concentration
    /// @param[out] totmob    total mobilities.
    void computeTotalMobility(const Opm::IncompPropertiesInterface& props,
			      const Opm::PolymerProperties& polyprops,
			      const std::vector<int>& cells,
			      const std::vector<double>& s,
			      const std::vector<double>& c,
			      const std::vector<double>& cmax,
			      std::vector<double>& totmob)
    {
	int num_cells = cells.size();
	totmob.resize(num_cells);
	std::vector<double> kr(2*num_cells);
	props.relperm(num_cells, &s[0], &cells[0], &kr[0], 0);
        const double* visc = props.viscosity();
	for (int cell = 0; cell < num_cells; ++cell) {
            double*  kr_cell = &kr[2*cell];
            polyprops.effectiveTotalMobility(c[cell], cmax[cell], visc, kr_cell,
                                             totmob[cell]);
	}
    }
    TransportSolverTwophaseReorder::TransportSolverTwophaseReorder(const UnstructuredGrid& grid,
                                                                   const Opm::IncompPropertiesInterface& props,
                                                                   const double* gravity,
                                                                   const double tol,
                                                                   const int maxit)
        : grid_(grid),
          props_(props),
          tol_(tol),
          maxit_(maxit),
          darcyflux_(0),
          source_(0),
          dt_(0.0),
          saturation_(grid.number_of_cells, -1.0),
          fractionalflow_(grid.number_of_cells, -1.0),
          reorder_iterations_(grid.number_of_cells, 0),
          mob_(2*grid.number_of_cells, -1.0)
#ifdef EXPERIMENT_GAUSS_SEIDEL
        , ia_upw_(grid.number_of_cells + 1, -1),
          ja_upw_(grid.number_of_faces, -1),
          ia_downw_(grid.number_of_cells + 1, -1),
          ja_downw_(grid.number_of_faces, -1)
#endif
    {
        if (props.numPhases() != 2) {
            OPM_THROW(std::runtime_error, "Property object must have 2 phases");
        }
        visc_ = props.viscosity();
        int num_cells = props.numCells();
        smin_.resize(props.numPhases()*num_cells);
        smax_.resize(props.numPhases()*num_cells);
        std::vector<int> cells(num_cells);
        for (int i = 0; i < num_cells; ++i) {
            cells[i] = i;
        }
        props.satRange(props.numCells(), &cells[0], &smin_[0], &smax_[0]);
        if (gravity) {
            initGravity(gravity);
            initColumns();
        }
    }
예제 #5
0
std::vector<ADB>
phaseMobility(const Opm::IncompPropertiesInterface& props,
              const std::vector<int>& cells,
              const typename ADB::V& sw)
{
    typedef Eigen::Array<double, Eigen::Dynamic, 2, Eigen::RowMajor> TwoCol;
    typedef Eigen::Array<double, Eigen::Dynamic, 4, Eigen::RowMajor> FourCol;
    typedef Eigen::SparseMatrix<double> S;
    typedef typename ADB::V V;
    typedef typename ADB::M M;
    const int nc = props.numCells();
    TwoCol s(nc, 2);
    s.leftCols<1>() = sw;
    s.rightCols<1>() = 1.0 - s.leftCols<1>();
    TwoCol kr(nc, 2);
    FourCol dkr(nc, 4);
    props.relperm(nc, s.data(), cells.data(), kr.data(), dkr.data());
    V krw = kr.leftCols<1>();
    V kro = kr.rightCols<1>();
    V dkrw = dkr.leftCols<1>();  // Left column is top-left of dkr/ds 2x2 matrix.
    V dkro = -dkr.rightCols<1>(); // Right column is bottom-right of dkr/ds 2x2 matrix.
    S krwjac(nc,nc);
    S krojac(nc,nc);
    auto sizes = Eigen::ArrayXi::Ones(nc);
    krwjac.reserve(sizes);
    krojac.reserve(sizes);
    for (int c = 0; c < nc; ++c) {
        krwjac.insert(c,c) = dkrw(c);
        krojac.insert(c,c) = dkro(c);
    }
    const double* mu = props.viscosity();
    std::vector<M> dmw = { M(krwjac)/mu[0] };
    std::vector<M> dmo = { M(krojac)/mu[1] };

    std::vector<ADB> pmobc = { ADB::function(krw / mu[0], std::move(dmw)) ,
                               ADB::function(kro / mu[1], std::move(dmo)) };
    return pmobc;
}
예제 #6
0
    /// @brief Computes phase mobilities for a set of saturation values.
    /// @param[in]  props     rock and fluid properties
    /// @param[in]  cells     cells with which the saturation values are associated
    /// @param[in]  s         saturation values (for all phases)
    /// @param[out] pmobc     phase mobilities (for all phases).
    void computePhaseMobilities(const Opm::IncompPropertiesInterface& props,
                                const std::vector<int>&               cells,
                                const std::vector<double>&            s    ,
                                std::vector<double>&                  pmobc)
    {
        const std::vector<int>::size_type nc = cells.size();
        const std::size_t                 np = props.numPhases();

        assert(s.size() == nc * np);

        std::vector<double>(nc * np, 0.0).swap(pmobc );
        double*                                dpmobc = 0;
        props.relperm(static_cast<const int>(nc), &s[0], &cells[0],
                      &pmobc[0], dpmobc);

        const double*                 mu  = props.viscosity();
        std::vector<double>::iterator lam = pmobc.begin();
        for (std::vector<int>::size_type c = 0; c < nc; ++c) {
            for (std::size_t p = 0; p < np; ++p, ++lam) {
                *lam /= mu[ p ];
            }
        }
    }