void TransportSolverTwophaseReorder::solve(const double* porevolume,
                                               const double* source,
                                               const double dt,
                                               TwophaseState& state)
    {
        darcyflux_ = &state.faceflux()[0];
        porevolume_ = porevolume;
        source_ = source;
        dt_ = dt;
        toWaterSat(state.saturation(), saturation_);

#ifdef EXPERIMENT_GAUSS_SEIDEL
        std::vector<int> seq(grid_.number_of_cells);
        std::vector<int> comp(grid_.number_of_cells + 1);
        int ncomp;
        compute_sequence_graph(&grid_, darcyflux_,
                               &seq[0], &comp[0], &ncomp,
                               &ia_upw_[0], &ja_upw_[0]);
        const int nf = grid_.number_of_faces;
        std::vector<double> neg_darcyflux(nf);
        std::transform(darcyflux_, darcyflux_ + nf, neg_darcyflux.begin(), std::negate<double>());
        compute_sequence_graph(&grid_, &neg_darcyflux[0],
                               &seq[0], &comp[0], &ncomp,
                               &ia_downw_[0], &ja_downw_[0]);
#endif
        std::fill(reorder_iterations_.begin(),reorder_iterations_.end(),0);
        reorderAndTransport(grid_, darcyflux_);
        toBothSat(saturation_, state.saturation());
    }
    void TransportSolverTwophaseReorder::solveGravity(const double* porevolume,
                                                      const double dt,
                                                      TwophaseState& state)
    {
        // Initialize mobilities.
        const int nc = grid_.number_of_cells;
        std::vector<int> cells(nc);
        for (int c = 0; c < nc; ++c) {
            cells[c] = c;
        }
        mob_.resize(2*nc);
        props_.relperm(cells.size(), &state.saturation()[0], &cells[0], &mob_[0], 0);
        const double* mu = props_.viscosity();
        for (int c = 0; c < nc; ++c) {
            mob_[2*c] /= mu[0];
            mob_[2*c + 1] /= mu[1];
        }

        // Set up other variables.
        porevolume_ = porevolume;
        dt_ = dt;
        toWaterSat(state.saturation(), saturation_);

        // Solve on all columns.
        int num_iters = 0;
        for (std::vector<std::vector<int> >::size_type i = 0; i < columns_.size(); i++) {
            // std::cout << "==== new column" << std::endl;
            num_iters += solveGravityColumn(columns_[i]);
        }
        std::cout << "Gauss-Seidel column solver average iterations: "
                  << double(num_iters)/double(columns_.size()) << std::endl;

        toBothSat(saturation_, state.saturation());
    }
Пример #3
0
void
VertEqImpl::downscale (const TwophaseState &coarseScale,
                       TwophaseState &fineScale) {
	// assume that the fineScale storage is already initialized
	if (!fineScale.pressure().size() == ts->number_of_cells) {
		throw OPM_EXC ("Fine scale state is not dimensioned correctly");
	}

	// properties object handle the actual downscaling since it
	// already has the information about the interface.
	// update the coarse saturation *before* we downscale to 3D,
	// since we need the residual interface for that.
	pr->upd_res_sat (&coarseScale.saturation ()[0]);
	pr->downscale_saturation (&coarseScale.saturation ()[0],
	                          &fineScale.saturation ()[0]);
	pr->downscale_pressure (&coarseScale.saturation ()[0],
	                        &coarseScale.pressure ()[0],
	                        &fineScale.pressure ()[0]);
}
Пример #4
0
    /// Compute per-solve dynamic properties.
    void IncompTpfa::computePerSolveDynamicData(const double /*dt*/,
                                                const TwophaseState& state,
                                                const WellState& /*well_state*/)
    {
        // Computed here:
        //
        // std::vector<double> wdp_;
        // std::vector<double> totmob_;
        // std::vector<double> omega_;
        // std::vector<double> trans_;
        // std::vector<double> gpress_omegaweighted_;
        // std::vector<double> initial_porevol_;
        // ifs_tpfa_forces forces_;

        // wdp_
        if (wells_) {
            Opm::computeWDP(*wells_, grid_, state.saturation(), props_.density(),
                            gravity_ ? gravity_[2] : 0.0, true, wdp_);
        }
        // totmob_, omega_, gpress_omegaweighted_
        if (gravity_) {
            computeTotalMobilityOmega(props_, allcells_, state.saturation(), totmob_, omega_);
            mim_ip_density_update(grid_.number_of_cells, grid_.cell_facepos,
                                  &omega_[0],
                                  &gpress_[0], &gpress_omegaweighted_[0]);
        } else {
            computeTotalMobility(props_, allcells_, state.saturation(), totmob_);
        }
        // trans_
        tpfa_eff_trans_compute(const_cast<UnstructuredGrid*>(&grid_), &totmob_[0], &htrans_[0], &trans_[0]);
        // initial_porevol_
        if (rock_comp_props_ && rock_comp_props_->isActive()) {
            computePorevolume(grid_, props_.porosity(), *rock_comp_props_, state.pressure(), initial_porevol_);
        }
        // forces_
        forces_.src = src_.empty() ? NULL : &src_[0];
        forces_.bc = bcs_;
        forces_.W = wells_;
        forces_.totmob = &totmob_[0];
        forces_.wdp = wdp_.empty() ? NULL : &wdp_[0];
    }
Пример #5
0
void
VertEqImpl::upscale (const TwophaseState& fineScale,
                     TwophaseState& coarseScale) {
	// dimension state object to the top grid
	coarseScale.init (*ts, pr->numPhases ());

	// upscale pressure and saturation to find the initial state of
	// the two-dimensional domain. we only need to set the pressure
	// and saturation, the flux is an output field. these methods
	// are handled by the props class, since it already has access to
	// the densities and weights.
	pr->upscale_saturation (&fineScale.saturation ()[0],
	                        &coarseScale.saturation ()[0]);
	pr->upd_res_sat (&coarseScale.saturation ()[0]);
	pr->upscale_pressure (&coarseScale.saturation ()[0],
	                      &fineScale.pressure ()[0],
	                      &coarseScale.pressure ()[0]);

	// use the regular helper method to initialize the face pressure
	// since it is implemented in the header, we have access to it
	// even though it is in an anonymous namespace!
	const UnstructuredGrid& g = this->grid();

	initFacePressure (UgGridHelpers::dimensions (g),
	                  UgGridHelpers::numFaces (g),
	                  UgGridHelpers::faceCells (g),
	                  UgGridHelpers::beginFaceCentroids (g),
	                  UgGridHelpers::beginCellCentroids (g),
	                  coarseScale);

	// update the properties from the initial state (the
	// simulation object won't call this method before the
	// first timestep; it assumes that the state is initialized
	// accordingly (which is what we do here now)
	notify (coarseScale);
}
Пример #6
0
/// \page tutorial3
/// \section commentedsource1 Program walk-through.
/// \details
/// Main function
/// \snippet tutorial3.cpp main
/// \internal [main]
int main ()
try
{
    /// \internal [main]
    /// \endinternal

    /// \page tutorial3
    /// \details
    /// We define the grid. A Cartesian grid with 400 cells,
    /// each being 10m along each side. Note that we treat the
    /// grid as 3-dimensional, but have a thickness of only one
    /// layer in the Z direction.
    ///
    /// The Opm::GridManager is responsible for creating and destroying the grid,
    /// the UnstructuredGrid data structure contains the actual grid topology
    /// and geometry.
    /// \snippet tutorial3.cpp grid
    /// \internal [grid]
    int nx = 20;
    int ny = 20;
    int nz = 1;
    double dx = 10.0;
    double dy = 10.0;
    double dz = 10.0;
    using namespace Opm;
    GridManager grid_manager(nx, ny, nz, dx, dy, dz);
    const UnstructuredGrid& grid = *grid_manager.c_grid();
    int num_cells = grid.number_of_cells;
    /// \internal [grid]
    /// \endinternal

    /// \page tutorial3
    /// \details
    /// We define the properties of the fluid.\n
    /// Number of phases, phase densities, phase viscosities,
    /// rock porosity and permeability.
    ///
    /// We always use SI units in the simulator. Many units are
    /// available for use, however.  They are stored as constants in
    /// the Opm::unit namespace, while prefixes are in the Opm::prefix
    /// namespace. See Units.hpp for more.
    /// \snippet tutorial3.cpp set properties
    /// \internal [set properties]
    int num_phases = 2;
    using namespace Opm::unit;
    using namespace Opm::prefix;
    std::vector<double> density(num_phases, 1000.0);
    std::vector<double> viscosity(num_phases, 1.0*centi*Poise);
    double porosity = 0.5;
    double permeability = 10.0*milli*darcy;
    /// \internal [set properties]
    /// \endinternal

    /// \page tutorial3
    /// \details We define the relative permeability function. We use a basic fluid
    /// description and set this function to be linear. For more realistic fluid, the
    /// saturation function may be interpolated from experimental data.
    /// \snippet tutorial3.cpp relperm
    /// \internal [relperm]
    SaturationPropsBasic::RelPermFunc rel_perm_func = SaturationPropsBasic::Linear;
    /// \internal [relperm]
    /// \endinternal

    /// \page tutorial3
    /// \details We construct a basic fluid and rock property object
    /// with the properties we have defined above.  Each property is
    /// constant and hold for all cells.
    /// \snippet tutorial3.cpp properties
    /// \internal [properties]
    IncompPropertiesBasic props(num_phases, rel_perm_func, density, viscosity,
                                porosity, permeability, grid.dimensions, num_cells);
    /// \internal [properties]
    /// \endinternal

    /// \page tutorial3
    /// \details Gravity parameters. Here, we set zero gravity.
    /// \snippet tutorial3.cpp gravity
    /// \internal [gravity]
    const double *grav = 0;
    std::vector<double> omega;
    /// \internal [gravity]
    /// \endinternal

    /// \page tutorial3
    /// \details We set up the source term. Positive numbers indicate that the cell is a source,
    /// while negative numbers indicate a sink.
    /// \snippet tutorial3.cpp source
    /// \internal [source]
    std::vector<double> src(num_cells, 0.0);
    src[0] = 1.;
    src[num_cells-1] = -1.;
    /// \internal [source]
    /// \endinternal

    /// \page tutorial3
    /// \details We set up the boundary conditions. Letting bcs be empty is equivalent
    /// to no-flow boundary conditions.
    /// \snippet tutorial3.cpp boundary
    /// \internal [boundary]
    FlowBCManager bcs;
    /// \internal [boundary]
    /// \endinternal

    /// \page tutorial3
    /// \details We may now set up the pressure solver. At this point,
    /// unchanging parameters such as transmissibility are computed
    /// and stored internally by the IncompTpfa class. The null pointer
    /// constructor argument is for wells, which are not used in this tutorial.
    /// \snippet tutorial3.cpp pressure solver
    /// \internal [pressure solver]
    LinearSolverUmfpack linsolver;
    IncompTpfa psolver(grid, props, linsolver, grav, NULL, src, bcs.c_bcs());
    /// \internal [pressure solver]
    /// \endinternal

    /// \page tutorial3
    /// \details We set up a state object for the wells. Here, there are
    /// no wells and we let it remain empty.
    /// \snippet tutorial3.cpp well
    /// \internal [well]
    WellState well_state;
    /// \internal [well]
    /// \endinternal

    /// \page tutorial3
    /// \details We compute the pore volume
    /// \snippet tutorial3.cpp pore volume
    /// \internal [pore volume]
    std::vector<double> porevol;
    Opm::computePorevolume(grid, props.porosity(), porevol);
    /// \internal [pore volume]
    /// \endinternal

    /// \page tutorial3
    /// \details Set up the transport solver. This is a reordering implicit Euler transport solver.
    /// \snippet tutorial3.cpp transport solver
    /// \internal [transport solver]
    const double tolerance = 1e-9;
    const int max_iterations = 30;
    Opm::TransportSolverTwophaseReorder transport_solver(grid, props, NULL, tolerance, max_iterations);
    /// \internal [transport solver]
    /// \endinternal

    /// \page tutorial3
    /// \details Time integration parameters
    /// \snippet tutorial3.cpp time parameters
    /// \internal [time parameters]
    const double dt = 0.1*day;
    const int num_time_steps = 20;
    /// \internal [time parameters]
    /// \endinternal


    /// \page tutorial3
    /// \details We define a vector which contains all cell indexes. We use this
    /// vector to set up parameters on the whole domain.
    /// \snippet tutorial3.cpp cell indexes
    /// \internal [cell indexes]
    std::vector<int> allcells(num_cells);
    for (int cell = 0; cell < num_cells; ++cell) {
        allcells[cell] = cell;
    }
    /// \internal [cell indexes]
    /// \endinternal

    /// \page tutorial3
    /// \details
    /// We set up a two-phase state object, and
    /// initialize water saturation to minimum everywhere.
    /// \snippet tutorial3.cpp two-phase state
    /// \internal [two-phase state]
    TwophaseState state;
    state.init(grid.number_of_cells , grid.number_of_faces, 2);
    initSaturation( allcells , props , state , MinSat );

    /// \internal [two-phase state]
    /// \endinternal

    /// \page tutorial3
    /// \details This string stream will be used to construct a new
    /// output filename at each timestep.
    /// \snippet tutorial3.cpp output stream
    /// \internal [output stream]
    std::ostringstream vtkfilename;
    /// \internal [output stream]
    /// \endinternal


    /// \page tutorial3
    /// \details Loop over the time steps.
    /// \snippet tutorial3.cpp time loop
    /// \internal [time loop]
    for (int i = 0; i < num_time_steps; ++i) {
        /// \internal [time loop]
	/// \endinternal


        /// \page tutorial3
        /// \details Solve the pressure equation
        /// \snippet tutorial3.cpp solve pressure
        /// \internal [solve pressure]
        psolver.solve(dt, state, well_state);
        /// \internal [solve pressure]
	/// \endinternal

        /// \page tutorial3
        /// \details  Solve the transport equation.
        /// \snippet tutorial3.cpp transport solve
	/// \internal [transport solve]
        transport_solver.solve(&porevol[0], &src[0], dt, state);
        /// \internal [transport solve]
	/// \endinternal

        /// \page tutorial3
        /// \details Write the output to file.
        /// \snippet tutorial3.cpp write output
	/// \internal [write output]
        vtkfilename.str("");
        vtkfilename << "tutorial3-" << std::setw(3) << std::setfill('0') << i << ".vtu";
        std::ofstream vtkfile(vtkfilename.str().c_str());
        Opm::DataMap dm;
        dm["saturation"] = &state.saturation();
        dm["pressure"] = &state.pressure();
        Opm::writeVtkData(grid, dm, vtkfile);
    }
}
catch (const std::exception &e) {
    std::cerr << "Program threw an exception: " << e.what() << "\n";
    throw;
}
Пример #7
0
void
VertEqImpl::notify (const TwophaseState& coarseScale) {
	// forward this request to the properties we have stored
	pr->upd_res_sat (&coarseScale.saturation()[0]);
}
Пример #8
0
    SimulatorReport SimulatorIncompTwophase::Impl::run(SimulatorTimer& timer,
                                                       TwophaseState& state,
                                                       WellState& well_state)
    {
        std::vector<double> transport_src;

        // Initialisation.
        std::vector<double> porevol;
        if (rock_comp_props_ && rock_comp_props_->isActive()) {
            computePorevolume(grid_, props_.porosity(), *rock_comp_props_, state.pressure(), porevol);
        } else {
            computePorevolume(grid_, props_.porosity(), porevol);
        }
        const double tot_porevol_init = std::accumulate(porevol.begin(), porevol.end(), 0.0);
        std::vector<double> initial_porevol = porevol;

        // Main simulation loop.
        Opm::time::StopWatch pressure_timer;
        double ptime = 0.0;
        Opm::time::StopWatch transport_timer;
        double ttime = 0.0;
        Opm::time::StopWatch callback_timer;
        double time_in_callbacks = 0.0;
        Opm::time::StopWatch step_timer;
        Opm::time::StopWatch total_timer;
        total_timer.start();
        double init_satvol[2] = { 0.0 };
        double satvol[2] = { 0.0 };
        double tot_injected[2] = { 0.0 };
        double tot_produced[2] = { 0.0 };
        Opm::computeSaturatedVol(porevol, state.saturation(), init_satvol);
        *log_ << "\nInitial saturations are    " << init_satvol[0]/tot_porevol_init
              << "    " << init_satvol[1]/tot_porevol_init << std::endl;
        Opm::Watercut watercut;
        watercut.push(0.0, 0.0, 0.0);
        Opm::WellReport wellreport;
        std::vector<double> fractional_flows;
        std::vector<double> well_resflows_phase;
        if (wells_) {
            well_resflows_phase.resize((wells_->number_of_phases)*(wells_->number_of_wells), 0.0);
            wellreport.push(props_, *wells_, state.saturation(), 0.0, well_state.bhp(), well_state.perfRates());
        }
        std::fstream tstep_os;
        if (output_) {
            std::string filename = output_dir_ + "/step_timing.param";
            tstep_os.open(filename.c_str(), std::fstream::out | std::fstream::app);
        }
        while (!timer.done()) {
            // Report timestep and (optionally) write state to disk.
            step_timer.start();
            timer.report(*log_);
            if (output_ && (timer.currentStepNum() % output_interval_ == 0)) {
                if (output_vtk_) {
                    outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
                }
                outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
                if (use_reorder_) {
                    // This use of dynamic_cast is not ideal, but should be safe.
                    outputVectorMatlab(std::string("reorder_it"),
                                       dynamic_cast<const TransportSolverTwophaseReorder&>(*tsolver_).getReorderIterations(),
                                       timer.currentStepNum(), output_dir_);
                }
            }

            SimulatorReport sreport;

            // Solve pressure equation.
            if (check_well_controls_) {
                computeFractionalFlow(props_, allcells_, state.saturation(), fractional_flows);
                wells_manager_.applyExplicitReinjectionControls(well_resflows_phase, well_resflows_phase);
            }
            bool well_control_passed = !check_well_controls_;
            int well_control_iteration = 0;
            do {
                // Run solver.
                pressure_timer.start();
                std::vector<double> initial_pressure = state.pressure();
                psolver_.solve(timer.currentStepLength(), state, well_state);

                // Renormalize pressure if rock is incompressible, and
                // there are no pressure conditions (bcs or wells).
                // It is deemed sufficient for now to renormalize
                // using geometric volume instead of pore volume.
                if ((rock_comp_props_ == NULL || !rock_comp_props_->isActive())
                    && allNeumannBCs(bcs_) && allRateWells(wells_)) {
                    // Compute average pressures of previous and last
                    // step, and total volume.
                    double av_prev_press = 0.0;
                    double av_press = 0.0;
                    double tot_vol = 0.0;
                    const int num_cells = grid_.number_of_cells;
                    for (int cell = 0; cell < num_cells; ++cell) {
                        av_prev_press += initial_pressure[cell]*grid_.cell_volumes[cell];
                        av_press      += state.pressure()[cell]*grid_.cell_volumes[cell];
                        tot_vol       += grid_.cell_volumes[cell];
                    }
                    // Renormalization constant
                    const double ren_const = (av_prev_press - av_press)/tot_vol;
                    for (int cell = 0; cell < num_cells; ++cell) {
                        state.pressure()[cell] += ren_const;
                    }
                    const int num_wells = (wells_ == NULL) ? 0 : wells_->number_of_wells;
                    for (int well = 0; well < num_wells; ++well) {
                        well_state.bhp()[well] += ren_const;
                    }
                }

                // Stop timer and report.
                pressure_timer.stop();
                double pt = pressure_timer.secsSinceStart();
                *log_ << "Pressure solver took:  " << pt << " seconds." << std::endl;
                ptime += pt;
                sreport.pressure_time = pt;

                // Optionally, check if well controls are satisfied.
                if (check_well_controls_) {
                    Opm::computePhaseFlowRatesPerWell(*wells_,
                                                      well_state.perfRates(),
                                                      fractional_flows,
                                                      well_resflows_phase);
                    *log_ << "Checking well conditions." << std::endl;
                    // For testing we set surface := reservoir
                    well_control_passed = wells_manager_.conditionsMet(well_state.bhp(), well_resflows_phase, well_resflows_phase);
                    ++well_control_iteration;
                    if (!well_control_passed && well_control_iteration > max_well_control_iterations_) {
                        OPM_THROW(std::runtime_error, "Could not satisfy well conditions in " << max_well_control_iterations_ << " tries.");
                    }
                    if (!well_control_passed) {
                        *log_ << "Well controls not passed, solving again." << std::endl;
                    } else {
                        *log_ << "Well conditions met." << std::endl;
                    }
                }
            } while (!well_control_passed);

            // Update pore volumes if rock is compressible.
            if (rock_comp_props_ && rock_comp_props_->isActive()) {
                initial_porevol = porevol;
                computePorevolume(grid_, props_.porosity(), *rock_comp_props_, state.pressure(), porevol);
            }

            // Process transport sources (to include bdy terms and well flows).
            Opm::computeTransportSource(grid_, src_, state.faceflux(), 1.0,
                                        wells_, well_state.perfRates(), transport_src);

            // Solve transport.
            transport_timer.start();
            double stepsize = timer.currentStepLength();
            if (num_transport_substeps_ != 1) {
                stepsize /= double(num_transport_substeps_);
                *log_ << "Making " << num_transport_substeps_ << " transport substeps." << std::endl;
            }
            double injected[2] = { 0.0 };
            double produced[2] = { 0.0 };
            for (int tr_substep = 0; tr_substep < num_transport_substeps_; ++tr_substep) {
                tsolver_->solve(&initial_porevol[0], &transport_src[0], stepsize, state);

                double substep_injected[2] = { 0.0 };
                double substep_produced[2] = { 0.0 };
                Opm::computeInjectedProduced(props_, state.saturation(), transport_src, stepsize,
                                             substep_injected, substep_produced);
                injected[0] += substep_injected[0];
                injected[1] += substep_injected[1];
                produced[0] += substep_produced[0];
                produced[1] += substep_produced[1];
                if (use_reorder_ && use_segregation_split_) {
                    // Again, unfortunate but safe use of dynamic_cast.
                    // Possible solution: refactor gravity solver to its own class.
                    dynamic_cast<TransportSolverTwophaseReorder&>(*tsolver_)
                        .solveGravity(&initial_porevol[0], stepsize, state);
                }
                watercut.push(timer.simulationTimeElapsed() + timer.currentStepLength(),
                              produced[0]/(produced[0] + produced[1]),
                              tot_produced[0]/tot_porevol_init);
                if (wells_) {
                    wellreport.push(props_, *wells_, state.saturation(),
                                    timer.simulationTimeElapsed() + timer.currentStepLength(),
                                    well_state.bhp(), well_state.perfRates());
                }
            }
            transport_timer.stop();
            double tt = transport_timer.secsSinceStart();
            sreport.transport_time = tt;
            *log_ << "Transport solver took: " << tt << " seconds." << std::endl;
            ttime += tt;
            // Report volume balances.
            Opm::computeSaturatedVol(porevol, state.saturation(), satvol);
            tot_injected[0] += injected[0];
            tot_injected[1] += injected[1];
            tot_produced[0] += produced[0];
            tot_produced[1] += produced[1];
            reportVolumes(*log_, satvol, tot_porevol_init,
                          tot_injected, tot_produced,
                          injected, produced,
                          init_satvol);
            sreport.total_time =  step_timer.secsSinceStart();
            if (output_) {
                sreport.reportParam(tstep_os);
            }

            // advance the timer to the end of the timestep *before* notifying
            // the client that the timestep is done
            ++timer;

            // notify all clients that we are done with the timestep
            callback_timer.start ();
            timestep_completed_.signal ();
            callback_timer.stop ();
            time_in_callbacks += callback_timer.secsSinceStart ();
        }

        if (output_) {
            if (output_vtk_) {
                outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
            }
            outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
            if (use_reorder_) {
                // This use of dynamic_cast is not ideal, but should be safe.
                outputVectorMatlab(std::string("reorder_it"),
                                   dynamic_cast<const TransportSolverTwophaseReorder&>(*tsolver_).getReorderIterations(),
                                   timer.currentStepNum(), output_dir_);
                }
            outputWaterCut(watercut, output_dir_);
            if (wells_) {
                outputWellReport(wellreport, output_dir_);
            }
            tstep_os.close();
        }

        total_timer.stop();

        SimulatorReport report;
        report.pressure_time = ptime;
        report.transport_time = ttime;
        report.total_time = total_timer.secsSinceStart() - time_in_callbacks;
        return report;
    }