/// @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;
            }
        }
    }
Esempio n. 2
0
    SimulatorReport SimulatorPolymer::Impl::run(SimulatorTimer& timer,
                                                PolymerState& 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);


        // Main simulation loop.
        Opm::time::StopWatch pressure_timer;
        double ptime = 0.0;
        Opm::time::StopWatch transport_timer;
        double ttime = 0.0;
        Opm::time::StopWatch total_timer;
        total_timer.start();
        double init_satvol[2] = { 0.0 };
        double init_polymass = 0.0;
        double satvol[2] = { 0.0 };
        double polymass = 0.0;
        double polymass_adsorbed = 0.0;
        double injected[2] = { 0.0 };
        double produced[2] = { 0.0 };
        double polyinj = 0.0;
        double polyprod = 0.0;
        double tot_injected[2] = { 0.0 };
        double tot_produced[2] = { 0.0 };
        double tot_polyinj = 0.0;
        double tot_polyprod = 0.0;
        Opm::computeSaturatedVol(porevol, state.saturation(), init_satvol);
        std::cout << "\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());
        }
        for (; !timer.done(); ++timer) {
            // Report timestep and (optionally) write state to disk.
            timer.report(std::cout);
            if (output_ && (timer.currentStepNum() % output_interval_ == 0)) {
                outputState(grid_, state, timer.currentStepNum(), output_dir_);
            }

            // Solve pressure.
            do {
                pressure_timer.start();
                psolver_.solve(timer.currentStepLength(), state, well_state);
                pressure_timer.stop();
                double pt = pressure_timer.secsSinceStart();
                std::cout << "Pressure solver took:  " << pt << " seconds." << std::endl;
                ptime += pt;
            } while (false);

            // Update pore volumes if rock is compressible.
            if (rock_comp_props_ && rock_comp_props_->isActive()) {
                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);

            // Find inflow rate.
            const double current_time = timer.currentTime();
            double stepsize = timer.currentStepLength();
            const double inflowc0 = poly_inflow_(current_time + 1e-5*stepsize);
            const double inflowc1 = poly_inflow_(current_time + (1.0 - 1e-5)*stepsize);
            if (inflowc0 != inflowc1) {
                std::cout << "**** Warning: polymer inflow rate changes during timestep. Using rate near start of step.";
            }
            const double inflow_c = inflowc0;

            // Solve transport.
            transport_timer.start();
            if (num_transport_substeps_ != 1) {
                stepsize /= double(num_transport_substeps_);
                std::cout << "Making " << num_transport_substeps_ << " transport substeps." << std::endl;
            }
            for (int tr_substep = 0; tr_substep < num_transport_substeps_; ++tr_substep) {
                tsolver_.solve(&state.faceflux()[0], &porevol[0], &transport_src[0], stepsize, inflow_c,
                               state.saturation(), state.concentration(), state.maxconcentration());
                Opm::computeInjectedProduced(props_, poly_props_,
                                             state.saturation(), state.concentration(), state.maxconcentration(),
                                             transport_src, timer.currentStepLength(), inflow_c,
                                             injected, produced, polyinj, polyprod);
                if (use_segregation_split_) {
                    tsolver_.solveGravity(columns_, &porevol[0], stepsize,
                                          state.saturation(), state.concentration(), state.maxconcentration());
                }
            }
            transport_timer.stop();
            double tt = transport_timer.secsSinceStart();
            std::cout << "Transport solver took: " << tt << " seconds." << std::endl;
            ttime += tt;

            // Report volume balances.
            Opm::computeSaturatedVol(porevol, state.saturation(), satvol);
            polymass = Opm::computePolymerMass(porevol, state.saturation(), state.concentration(), poly_props_.deadPoreVol());
            polymass_adsorbed = Opm::computePolymerAdsorbed(props_, poly_props_, porevol, state.maxconcentration());
            tot_injected[0] += injected[0];
            tot_injected[1] += injected[1];
            tot_produced[0] += produced[0];
            tot_produced[1] += produced[1];
            tot_polyinj += polyinj;
            tot_polyprod += polyprod;
            std::cout.precision(5);
            const int width = 18;
            std::cout << "\nVolume and polymer mass balance: "
                "   water(pv)           oil(pv)       polymer(kg)\n";
            std::cout << "    Saturated volumes:     "
                      << std::setw(width) << satvol[0]/tot_porevol_init
                      << std::setw(width) << satvol[1]/tot_porevol_init
                      << std::setw(width) << polymass << std::endl;
            std::cout << "    Adsorbed volumes:      "
                      << std::setw(width) << 0.0
                      << std::setw(width) << 0.0
                      << std::setw(width) << polymass_adsorbed << std::endl;
            std::cout << "    Injected volumes:      "
                      << std::setw(width) << injected[0]/tot_porevol_init
                      << std::setw(width) << injected[1]/tot_porevol_init
                      << std::setw(width) << polyinj << std::endl;
            std::cout << "    Produced volumes:      "
                      << std::setw(width) << produced[0]/tot_porevol_init
                      << std::setw(width) << produced[1]/tot_porevol_init
                      << std::setw(width) << polyprod << std::endl;
            std::cout << "    Total inj volumes:     "
                      << std::setw(width) << tot_injected[0]/tot_porevol_init
                      << std::setw(width) << tot_injected[1]/tot_porevol_init
                      << std::setw(width) << tot_polyinj << std::endl;
            std::cout << "    Total prod volumes:    "
                      << std::setw(width) << tot_produced[0]/tot_porevol_init
                      << std::setw(width) << tot_produced[1]/tot_porevol_init
                      << std::setw(width) << tot_polyprod << std::endl;
            std::cout << "    In-place + prod - inj: "
                      << std::setw(width) << (satvol[0] + tot_produced[0] - tot_injected[0])/tot_porevol_init
                      << std::setw(width) << (satvol[1] + tot_produced[1] - tot_injected[1])/tot_porevol_init
                      << std::setw(width) << (polymass + tot_polyprod - tot_polyinj + polymass_adsorbed) << std::endl;
            std::cout << "    Init - now - pr + inj: "
                      << std::setw(width) << (init_satvol[0] - satvol[0] - tot_produced[0] + tot_injected[0])/tot_porevol_init
                      << std::setw(width) << (init_satvol[1] - satvol[1] - tot_produced[1] + tot_injected[1])/tot_porevol_init
                      << std::setw(width) << (init_polymass - polymass - tot_polyprod + tot_polyinj - polymass_adsorbed)
                      << std::endl;
            std::cout.precision(8);

            watercut.push(timer.currentTime() + timer.currentStepLength(),
                          produced[0]/(produced[0] + produced[1]),
                          tot_produced[0]/tot_porevol_init);
            if (wells_) {
                wellreport.push(props_, *wells_, state.saturation(),
                                timer.currentTime() + timer.currentStepLength(),
                                well_state.bhp(), well_state.perfRates());
            }
        }

        if (output_) {
            outputState(grid_, state, timer.currentStepNum(), output_dir_);
            outputWaterCut(watercut, output_dir_);
            if (wells_) {
                outputWellReport(wellreport, output_dir_);
            }
        }

        total_timer.stop();

        SimulatorReport report;
        report.pressure_time = ptime;
        report.transport_time = ttime;
        report.total_time = total_timer.secsSinceStart();
        return report;
    }
// ----------------- Main program -----------------
int
main(int argc, char** argv)
try
{
    using namespace Opm;

    std::cout << "\n================    Test program for incompressible two-phase flow with polymer    ===============\n\n";
    parameter::ParameterGroup param(argc, argv, false);
    std::cout << "---------------    Reading parameters     ---------------" << std::endl;

    // If we have a "deck_filename", grid and props will be read from that.
    bool use_deck = param.has("deck_filename");
    boost::scoped_ptr<EclipseGridParser> deck;
    boost::scoped_ptr<GridManager> grid;
    boost::scoped_ptr<IncompPropertiesInterface> props;
    boost::scoped_ptr<RockCompressibility> rock_comp;
    PolymerState state;
    Opm::PolymerProperties poly_props;
    // bool check_well_controls = false;
    // int max_well_control_iterations = 0;
    double gravity[3] = { 0.0 };
    if (use_deck) {
        std::string deck_filename = param.get<std::string>("deck_filename");
        deck.reset(new EclipseGridParser(deck_filename));
        // Grid init
        grid.reset(new GridManager(*deck));
        // Rock and fluid init
        props.reset(new IncompPropertiesFromDeck(*deck, *grid->c_grid()));
        // check_well_controls = param.getDefault("check_well_controls", false);
        // max_well_control_iterations = param.getDefault("max_well_control_iterations", 10);
        // Rock compressibility.
        rock_comp.reset(new RockCompressibility(*deck));
        // Gravity.
        gravity[2] = deck->hasField("NOGRAV") ? 0.0 : unit::gravity;
        // Init state variables (saturation and pressure).
        if (param.has("init_saturation")) {
            initStateBasic(*grid->c_grid(), *props, param, gravity[2], state);
        } else {
            initStateFromDeck(*grid->c_grid(), *props, *deck, gravity[2], state);
        }
        // Init polymer properties.
        poly_props.readFromDeck(*deck);
    } else {
        // Grid init.
        const int nx = param.getDefault("nx", 100);
        const int ny = param.getDefault("ny", 100);
        const int nz = param.getDefault("nz", 1);
        const double dx = param.getDefault("dx", 1.0);
        const double dy = param.getDefault("dy", 1.0);
        const double dz = param.getDefault("dz", 1.0);
        grid.reset(new GridManager(nx, ny, nz, dx, dy, dz));
        // Rock and fluid init.
        props.reset(new IncompPropertiesBasic(param, grid->c_grid()->dimensions, grid->c_grid()->number_of_cells));
        // Rock compressibility.
        rock_comp.reset(new RockCompressibility(param));
        // Gravity.
        gravity[2] = param.getDefault("gravity", 0.0);
        // Init state variables (saturation and pressure).
        initStateBasic(*grid->c_grid(), *props, param, gravity[2], state);
        // Init Polymer state
        if (param.has("poly_init")) {
            double poly_init = param.getDefault("poly_init", 0.0);
            for (int cell = 0; cell < grid->c_grid()->number_of_cells; ++cell) {
                double smin[2], smax[2];
                props->satRange(1, &cell, smin, smax);
                if (state.saturation()[2*cell] > 0.5*(smin[0] + smax[0])) {
                    state.concentration()[cell] = poly_init;
                    state.maxconcentration()[cell] = poly_init;
                } else {
                    state.saturation()[2*cell + 0] = 0.;
                    state.saturation()[2*cell + 1] = 1.;
                    state.concentration()[cell] = 0.;
                    state.maxconcentration()[cell] = 0.;
                }
            }
        }
        // Init polymer properties.
        // Setting defaults to provide a simple example case.
        double c_max = param.getDefault("c_max_limit", 5.0);
        double mix_param = param.getDefault("mix_param", 1.0);
        double rock_density = param.getDefault("rock_density", 1000.0);
        double dead_pore_vol = param.getDefault("dead_pore_vol", 0.15);
        double res_factor = param.getDefault("res_factor", 1.) ; // res_factor = 1 gives no change in permeability
        double c_max_ads = param.getDefault("c_max_ads", 1.);
        int ads_index = param.getDefault<int>("ads_index", Opm::PolymerProperties::NoDesorption);
        std::vector<double> c_vals_visc(2, -1e100);
        c_vals_visc[0] = 0.0;
        c_vals_visc[1] = 7.0;
        std::vector<double> visc_mult_vals(2, -1e100);
        visc_mult_vals[0] = 1.0;
        // poly_props.visc_mult_vals[1] = param.getDefault("c_max_viscmult", 30.0);
        visc_mult_vals[1] = 20.0;
        std::vector<double> c_vals_ads(3, -1e100);
        c_vals_ads[0] = 0.0;
        c_vals_ads[1] = 2.0;
        c_vals_ads[2] = 8.0;
        std::vector<double> ads_vals(3, -1e100);
        ads_vals[0] = 0.0;
        ads_vals[1] = 0.0015;
        ads_vals[2] = 0.0025;
        // ads_vals[1] = 0.0;
        // ads_vals[2] = 0.0;
        poly_props.set(c_max, mix_param, rock_density, dead_pore_vol, res_factor, c_max_ads,
                       static_cast<Opm::PolymerProperties::AdsorptionBehaviour>(ads_index),
                       c_vals_visc,  visc_mult_vals, c_vals_ads, ads_vals);
    }

    // Warn if gravity but no density difference.
    bool use_gravity = (gravity[0] != 0.0 || gravity[1] != 0.0 || gravity[2] != 0.0);
    if (use_gravity) {
        if (props->density()[0] == props->density()[1]) {
            std::cout << "**** Warning: nonzero gravity, but zero density difference." << std::endl;
        }
    }
    const double *grav = use_gravity ? &gravity[0] : 0;

    // Initialising src
    int num_cells = grid->c_grid()->number_of_cells;
    std::vector<double> src(num_cells, 0.0);
    if (use_deck) {
        // Do nothing, wells will be the driving force, not source terms.
    } else {
        // Compute pore volumes, in order to enable specifying injection rate
        // terms of total pore volume.
        std::vector<double> porevol;
        if (rock_comp->isActive()) {
            computePorevolume(*grid->c_grid(), props->porosity(), *rock_comp, state.pressure(), porevol);
        } else {
            computePorevolume(*grid->c_grid(), props->porosity(), porevol);
        }
        const double tot_porevol_init = std::accumulate(porevol.begin(), porevol.end(), 0.0);
        const double default_injection = use_gravity ? 0.0 : 0.1;
        const double flow_per_sec = param.getDefault<double>("injected_porevolumes_per_day", default_injection)
            *tot_porevol_init/unit::day;
        src[0] = flow_per_sec;
        src[num_cells - 1] = -flow_per_sec;
    }

    // Boundary conditions.
    FlowBCManager bcs;
    if (param.getDefault("use_pside", false)) {
        int pside = param.get<int>("pside");
        double pside_pressure = param.get<double>("pside_pressure");
        bcs.pressureSide(*grid->c_grid(), FlowBCManager::Side(pside), pside_pressure);
    }

    // Linear solver.
    LinearSolverFactory linsolver(param);

    // Write parameters used for later reference.
    bool output = param.getDefault("output", true);
    if (output) {
      std::string output_dir =
        param.getDefault("output_dir", std::string("output"));
      boost::filesystem::path fpath(output_dir);
      try {
        create_directories(fpath);
      }
      catch (...) {
        OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath);
      }
      param.writeParam(output_dir + "/simulation.param");
    }


    std::cout << "\n\n================    Starting main simulation loop     ===============\n"
              << "                        (number of epochs: "
              << (use_deck ? deck->numberOfEpochs() : 1) << ")\n\n" << std::flush;

    SimulatorReport rep;
    if (!use_deck) {
        // Simple simulation without a deck.
        PolymerInflowBasic polymer_inflow(param.getDefault("poly_start_days", 300.0)*Opm::unit::day,
                                          param.getDefault("poly_end_days", 800.0)*Opm::unit::day,
                                          param.getDefault("poly_amount", poly_props.cMax()));
        WellsManager wells;
        SimulatorPolymer simulator(param,
                                   *grid->c_grid(),
                                   *props,
                                   poly_props,
                                   rock_comp->isActive() ? rock_comp.get() : 0,
                                   wells,
                                   polymer_inflow,
                                   src,
                                   bcs.c_bcs(),
                                   linsolver,
                                   grav);
        SimulatorTimer simtimer;
        simtimer.init(param);
        warnIfUnusedParams(param);
        WellState well_state;
        well_state.init(0, state);
        rep = simulator.run(simtimer, state, well_state);
    } else {
        // With a deck, we may have more epochs etc.
        WellState well_state;
        int step = 0;
        SimulatorTimer simtimer;
        // Use timer for last epoch to obtain total time.
        deck->setCurrentEpoch(deck->numberOfEpochs() - 1);
        simtimer.init(*deck);
        const double total_time = simtimer.totalTime();
        // Check for WPOLYMER presence in last epoch to decide
        // polymer injection control type.
        const bool use_wpolymer = deck->hasField("WPOLYMER");
        if (use_wpolymer) {
            if (param.has("poly_start_days")) {
                OPM_MESSAGE("Warning: Using WPOLYMER to control injection since it was found in deck. "
                        "You seem to be trying to control it via parameter poly_start_days (etc.) as well.");
            }
        }
        for (int epoch = 0; epoch < deck->numberOfEpochs(); ++epoch) {
            // Set epoch index.
            deck->setCurrentEpoch(epoch);

            // Update the timer.
            if (deck->hasField("TSTEP")) {
                simtimer.init(*deck);
            } else {
                if (epoch != 0) {
                    OPM_THROW(std::runtime_error, "No TSTEP in deck for epoch " << epoch);
                }
                simtimer.init(param);
            }
            simtimer.setCurrentStepNum(step);
            simtimer.setTotalTime(total_time);

            // Report on start of epoch.
            std::cout << "\n\n--------------    Starting epoch " << epoch << "    --------------"
                      << "\n                  (number of steps: "
                      << simtimer.numSteps() - step << ")\n\n" << std::flush;

            // Create new wells, polymer inflow controls.
            WellsManager wells(*deck, *grid->c_grid(), props->permeability());
            boost::scoped_ptr<PolymerInflowInterface> polymer_inflow;
            if (use_wpolymer) {
                if (wells.c_wells() == 0) {
                    OPM_THROW(std::runtime_error, "Cannot control polymer injection via WPOLYMER without wells.");
                }
                polymer_inflow.reset(new PolymerInflowFromDeck(*deck, *wells.c_wells(), props->numCells()));
            } else {
                polymer_inflow.reset(new PolymerInflowBasic(param.getDefault("poly_start_days", 300.0)*Opm::unit::day,
                                                            param.getDefault("poly_end_days", 800.0)*Opm::unit::day,
                                                            param.getDefault("poly_amount", poly_props.cMax())));
            }

            // @@@ HACK: we should really make a new well state and
            // properly transfer old well state to it every epoch,
            // since number of wells may change etc.
            if (epoch == 0) {
                well_state.init(wells.c_wells(), state);
            }

            // Create and run simulator.
            SimulatorPolymer simulator(param,
                                       *grid->c_grid(),
                                       *props,
                                       poly_props,
                                       rock_comp->isActive() ? rock_comp.get() : 0,
                                       wells,
                                       *polymer_inflow,
                                       src,
                                       bcs.c_bcs(),
                                       linsolver,
                                       grav);
            if (epoch == 0) {
                warnIfUnusedParams(param);
            }
            SimulatorReport epoch_rep = simulator.run(simtimer, state, well_state);

            // Update total timing report and remember step number.
            rep += epoch_rep;
            step = simtimer.currentStepNum();
        }
    }

    std::cout << "\n\n================    End of simulation     ===============\n\n";
    rep.report(std::cout);
}
catch (const std::exception &e) {
    std::cerr << "Program threw an exception: " << e.what() << "\n";
    throw;
}
Esempio n. 4
0
    SimulatorReport SimulatorPolymer::Impl::run(SimulatorTimer& timer,
                                                PolymerState& state,
                                                WellState& well_state)
    {
        std::vector<double> transport_src(grid_.number_of_cells);
        std::vector<double> polymer_inflow_c(grid_.number_of_cells);

        // 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 total_timer;
        total_timer.start();
        double init_satvol[2] = { 0.0 };
        double satvol[2] = { 0.0 };
        double polymass = computePolymerMass(porevol, state.saturation(), state.getCellData( state.CONCENTRATION ), poly_props_.deadPoreVol());
        double polymass_adsorbed = computePolymerAdsorbed(props_, poly_props_, porevol, state.getCellData( state.CMAX ));
        double init_polymass = polymass + polymass_adsorbed;
        double injected[2] = { 0.0 };
        double produced[2] = { 0.0 };
        double polyinj = 0.0;
        double polyprod = 0.0;
        double tot_injected[2] = { 0.0 };
        double tot_produced[2] = { 0.0 };
        double tot_polyinj = 0.0;
        double tot_polyprod = 0.0;
        Opm::computeSaturatedVol(porevol, state.saturation(), init_satvol);
        std::cout << "\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());
        }
        // Report timestep and (optionally) write state to disk.
        timer.report(std::cout);
        if (output_ && (timer.currentStepNum() % output_interval_ == 0)) {
            if (output_vtk_) {
                outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
            }
            if (output_binary_) {
                outputStateBinary(grid_, state, timer, output_dir_);
            }
            outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
        }

        // Solve pressure.
        if (check_well_controls_) {
            computeFractionalFlow(props_, poly_props_, allcells_,
                                  state.saturation(), state.getCellData( state.CONCENTRATION ), state.getCellData( state.CMAX ),
                                  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();
            std::cout << "Pressure solver took:  " << pt << " seconds." << std::endl;
            ptime += pt;

            // Optionally, check if well controls are satisfied.
            if (check_well_controls_) {
                Opm::computePhaseFlowRatesPerWell(*wells_,
                                                  well_state.perfRates(),
                                                  fractional_flows,
                                                  well_resflows_phase);
                std::cout << "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) {
                    std::cout << "Well controls not passed, solving again." << std::endl;
                } else {
                    std::cout << "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);

        // Find inflow rate.
        const double current_time = timer.simulationTimeElapsed();
        double stepsize = timer.currentStepLength();
        polymer_inflow_.getInflowValues(current_time, current_time + stepsize, polymer_inflow_c);

        // Solve transport.
        transport_timer.start();
        if (num_transport_substeps_ != 1) {
            stepsize /= double(num_transport_substeps_);
            std::cout << "Making " << num_transport_substeps_ << " transport substeps." << std::endl;
        }
        double substep_injected[2] = { 0.0 };
        double substep_produced[2] = { 0.0 };
        double substep_polyinj = 0.0;
        double substep_polyprod = 0.0;
        injected[0] = injected[1] = produced[0] = produced[1] = polyinj = polyprod = 0.0;
        for (int tr_substep = 0; tr_substep < num_transport_substeps_; ++tr_substep) {
            tsolver_.solve(&state.faceflux()[0], &initial_porevol[0], &transport_src[0], &polymer_inflow_c[0], stepsize,
                           state.saturation(), state.getCellData( state.CONCENTRATION ), state.getCellData( state.CMAX ));
            Opm::computeInjectedProduced(props_, poly_props_,
                                         state,
                                         transport_src, polymer_inflow_c, stepsize,
                                         substep_injected, substep_produced, substep_polyinj, substep_polyprod);
            injected[0] += substep_injected[0];
            injected[1] += substep_injected[1];
            produced[0] += substep_produced[0];
            produced[1] += substep_produced[1];
            polyinj += substep_polyinj;
            polyprod += substep_polyprod;
            if (use_segregation_split_) {
                tsolver_.solveGravity(columns_, &porevol[0], stepsize,
                                      state.saturation(), state.getCellData( state.CONCENTRATION ), state.getCellData( state.CMAX ));
            }
        }
        transport_timer.stop();
        double tt = transport_timer.secsSinceStart();
        std::cout << "Transport solver took: " << tt << " seconds." << std::endl;
        ttime += tt;

        // Report volume balances.
        Opm::computeSaturatedVol(porevol, state.saturation(), satvol);
        polymass = Opm::computePolymerMass(porevol, state.saturation(), state.getCellData( state.CONCENTRATION ), poly_props_.deadPoreVol());
        polymass_adsorbed = Opm::computePolymerAdsorbed(props_, poly_props_, porevol, state.getCellData( state.CMAX ));
        tot_injected[0] += injected[0];
        tot_injected[1] += injected[1];
        tot_produced[0] += produced[0];
        tot_produced[1] += produced[1];
        tot_polyinj += polyinj;
        tot_polyprod += polyprod;
        std::cout.precision(5);
        const int width = 18;
        std::cout << "\nVolume and polymer mass balance: "
            "   water(pv)           oil(pv)       polymer(kg)\n";
        std::cout << "    Saturated volumes:     "
                  << std::setw(width) << satvol[0]/tot_porevol_init
                  << std::setw(width) << satvol[1]/tot_porevol_init
                  << std::setw(width) << polymass << std::endl;
        std::cout << "    Adsorbed volumes:      "
                  << std::setw(width) << 0.0
                  << std::setw(width) << 0.0
                  << std::setw(width) << polymass_adsorbed << std::endl;
        std::cout << "    Injected volumes:      "
                  << std::setw(width) << injected[0]/tot_porevol_init
                  << std::setw(width) << injected[1]/tot_porevol_init
                  << std::setw(width) << polyinj << std::endl;
        std::cout << "    Produced volumes:      "
                  << std::setw(width) << produced[0]/tot_porevol_init
                  << std::setw(width) << produced[1]/tot_porevol_init
                  << std::setw(width) << polyprod << std::endl;
        std::cout << "    Total inj volumes:     "
                  << std::setw(width) << tot_injected[0]/tot_porevol_init
                  << std::setw(width) << tot_injected[1]/tot_porevol_init
                  << std::setw(width) << tot_polyinj << std::endl;
        std::cout << "    Total prod volumes:    "
                  << std::setw(width) << tot_produced[0]/tot_porevol_init
                  << std::setw(width) << tot_produced[1]/tot_porevol_init
                  << std::setw(width) << tot_polyprod << std::endl;
        std::cout << "    In-place + prod - inj: "
                  << std::setw(width) << (satvol[0] + tot_produced[0] - tot_injected[0])/tot_porevol_init
                  << std::setw(width) << (satvol[1] + tot_produced[1] - tot_injected[1])/tot_porevol_init
                  << std::setw(width) << (polymass + tot_polyprod - tot_polyinj + polymass_adsorbed) << std::endl;
        std::cout << "    Init - now - pr + inj: "
                  << std::setw(width) << (init_satvol[0] - satvol[0] - tot_produced[0] + tot_injected[0])/tot_porevol_init
                  << std::setw(width) << (init_satvol[1] - satvol[1] - tot_produced[1] + tot_injected[1])/tot_porevol_init
                  << std::setw(width) << (init_polymass - polymass - tot_polyprod + tot_polyinj - polymass_adsorbed)
                  << std::endl;
        std::cout.precision(8);

        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());
        }

        if (output_) {
            if (output_vtk_) {
                outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
            }
            if (output_binary_) {
                outputStateBinary(grid_, state, timer, output_dir_);
            }
            outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
            outputWaterCut(watercut, output_dir_);
            if (wells_) {
                outputWellReport(wellreport, output_dir_);
            }
        }

        total_timer.stop();

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