/* Copyright 2012 SINTEF ICT, Applied Mathematics. Copyright 2012 Statoil ASA. This file is part of the Open Porous Media Project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #if HAVE_CONFIG_H #include "config.h" #endif // HAVE_CONFIG_H #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include // #define EXPERIMENT_ISTL #ifdef EXPERIMENT_ISTL #include #endif #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include class AdHocProps : public Opm::IncompPropertiesBasic { public: AdHocProps(const Opm::parameter::ParameterGroup& param, int dim, int num_cells) : Opm::IncompPropertiesBasic(param, dim, num_cells) { ASSERT(numPhases() == 2); sw_.resize(3); sw_[0] = 0.2; sw_[1] = 0.7; sw_[2] = 1.0; krw_.resize(3); krw_[0] = 0.0; krw_[1] = 0.7; krw_[2] = 1.0; so_.resize(2); so_[0] = 0.3; so_[1] = 0.8; kro_.resize(2); kro_[0] = 0.0; kro_[1] = 1.0; } virtual void relperm(const int n, const double* s, const int* /*cells*/, double* kr, double* dkrds) const { // ASSERT(dkrds == 0); // We assume two phases flow for (int i = 0; i < n; ++i) { kr[2*i] = krw(s[2*i]); kr[2*i+1] = kro(s[2*i+1]); if (dkrds != 0) { dkrds[2*i] = krw_dsw(s[2*i]); dkrds[2*i+3] = kro_dso(s[2*i+1]); dkrds[2*i+1] = -dkrds[2*i+3]; dkrds[2*i+2] = -dkrds[2*i]; } } } virtual void satRange(const int n, const int* /*cells*/, double* smin, double* smax) const { const int np = 2; for (int i = 0; i < n; ++i) { smin[np*i + 0] = sw_[0]; smax[np*i + 0] = sw_.back(); smin[np*i + 1] = 1.0 - sw_[0]; smax[np*i + 1] = 1.0 - sw_.back(); } } private: double krw(double s) const { return Opm::linearInterpolation(sw_, krw_, s); } double krw_dsw(double s) const { return Opm::linearInterpolationDerivative(sw_, krw_, s); } double kro(double s) const { return Opm::linearInterpolation(so_, kro_, s); } double kro_dso(double s) const { return Opm::linearInterpolationDerivative(so_, kro_, s); } std::vector sw_; std::vector krw_; std::vector so_; std::vector kro_; }; // --------------- Types needed to define transport solver --------------- class PolymerFluid2pWrappingProps { public: PolymerFluid2pWrappingProps(const Opm::IncompPropertiesInterface& props, const Opm::PolymerProperties& polyprops) : props_(props), polyprops_(polyprops), smin_(props.numCells()*props.numPhases()), smax_(props.numCells()*props.numPhases()) { if (props.numPhases() != 2) { THROW("PolymerFluid2pWrapper requires 2 phases."); } const int num_cells = props.numCells(); std::vector cells(num_cells); for (int c = 0; c < num_cells; ++c) { cells[c] = c; } props.satRange(num_cells, &cells[0], &smin_[0], &smax_[0]); } double density(int phase) const { return props_.density()[phase]; } template void adsorption(const PolyC& c, const PolyC& cmax, CAds& cads, DCAdsDc& dcadsdc) { cads = polyprops_.adsorptionWithDer(c, cmax, &dcadsdc); } const std::vector porosity() const { const int num_cells = props_.numCells(); std::vector porosity(num_cells, 0.); const double* poro = props_.porosity(); for (std::vector::iterator it = porosity.begin(); it != porosity.end(); ++it, ++poro) { *it = poro[0]; } return porosity; } double rockdensity() const { return polyprops_.rockDensity(); } template void mobility(int cell, const Sat& s, const PolyC& c, Mob& mob, DMobDs& dmobds, DMobWatDc& dmobwatdc) const { const double* visc = props_.viscosity(); double relperm[2]; double drelpermds[4]; props_.relperm(1, &s[0], &cell, relperm, drelpermds); polyprops_.effectiveMobilities(c, visc, relperm, drelpermds, mob, dmobds, dmobwatdc); } template void pc(int c, const Sat& s, Pcap& pcap, DPcap& dpcap) const { double pc[2]; double dpc[4]; props_.capPress(1, &s[0], &c, pc, dpc); pcap = pc[0]; ASSERT(pc[1] == 0.0); dpcap = dpc[0]; ASSERT(dpc[1] == 0.0); ASSERT(dpc[2] == 0.0); ASSERT(dpc[3] == 0.0); } double s_min(int c) const { return smin_[2*c + 0]; } double s_max(int c) const { return smax_[2*c + 0]; } template void mc(const PolyC& c, Mc& mc, DMcDc& dmcdc) const { polyprops_.computeMc(c, mc, dmcdc); } private: const Opm::IncompPropertiesInterface& props_; const Opm::PolymerProperties& polyprops_; std::vector smin_; std::vector smax_; }; typedef PolymerFluid2pWrappingProps TwophaseFluidPolymer; typedef Opm::SinglePointUpwindTwoPhasePolymer NewtonPolymerTransportModel; using namespace Opm::ImplicitTransportDefault; typedef NewtonVectorCollection< ::std::vector > NVecColl; typedef JacobianSystem < struct CSRMatrix, NVecColl > JacSys; template class MaxNorm { public: static double norm(const Vector& v) { return AccumulationNorm ::norm(v); } }; typedef Opm::ImplicitTransport TransportSolver; static void outputState(const UnstructuredGrid& grid, const Opm::PolymerState& state, const int step, const std::string& output_dir) { // Write data in VTK format. std::ostringstream vtkfilename; vtkfilename << output_dir << "/output-" << std::setw(3) << std::setfill('0') << step << ".vtu"; std::ofstream vtkfile(vtkfilename.str().c_str()); if (!vtkfile) { THROW("Failed to open " << vtkfilename.str()); } Opm::DataMap dm; dm["saturation"] = &state.saturation(); dm["pressure"] = &state.pressure(); dm["concentration"] = &state.concentration(); dm["cmax"] = &state.maxconcentration(); std::vector cell_velocity; Opm::estimateCellVelocity(grid, state.faceflux(), cell_velocity); dm["velocity"] = &cell_velocity; Opm::writeVtkData(grid, dm, vtkfile); // Write data (not grid) in Matlab format for (Opm::DataMap::const_iterator it = dm.begin(); it != dm.end(); ++it) { std::ostringstream fname; fname << output_dir << "/" << it->first << "-" << std::setw(3) << std::setfill('0') << step << ".dat"; std::ofstream file(fname.str().c_str()); if (!file) { THROW("Failed to open " << fname.str()); } const std::vector& d = *(it->second); std::copy(d.begin(), d.end(), std::ostream_iterator(file, "\n")); } } class PolymerInflow { public: PolymerInflow(const double starttime, const double endtime, const double amount) : stime_(starttime), etime_(endtime), amount_(amount) { } double operator()(double time) { if (time >= stime_ && time < etime_) { return amount_; } else { return 0.0; } } private: double stime_; double etime_; double amount_; }; static void outputWaterCut(const Opm::Watercut& watercut, const std::string& output_dir) { // Write water cut curve. std::string fname = output_dir + "/watercut.txt"; std::ofstream os(fname.c_str()); if (!os) { THROW("Failed to open " << fname); } watercut.write(os); } // ----------------- Main program ----------------- int main(int argc, char** argv) { std::cout << "\n================ Test program for incompressible two-phase flow with polymer ===============\n\n"; Opm::parameter::ParameterGroup param(argc, argv, false); std::cout << "--------------- Reading parameters ---------------" << std::endl; // Reading various control parameters. const bool guess_old_solution = param.getDefault("guess_old_solution", false); const bool use_reorder = param.getDefault("use_reorder", true); const bool output = param.getDefault("output", true); std::string output_dir; int output_interval = 1; if (output) { output_dir = param.getDefault("output_dir", std::string("output")); // Ensure that output dir exists boost::filesystem::path fpath(output_dir); create_directories(fpath); output_interval = param.getDefault("output_interval", output_interval); } // If we have a "deck_filename", grid and props will be read from that. bool use_deck = param.has("deck_filename"); boost::scoped_ptr grid; boost::scoped_ptr props; boost::scoped_ptr wells; Opm::SimulatorTimer simtimer; Opm::PolymerState state; Opm::PolymerProperties polyprop; double gravity[3] = { 0.0 }; if (use_deck) { std::string deck_filename = param.get("deck_filename"); Opm::EclipseGridParser deck(deck_filename); // Grid init grid.reset(new Opm::GridManager(deck)); // Rock and fluid init const int* gc = grid->c_grid()->global_cell; std::vector global_cell(gc, gc + grid->c_grid()->number_of_cells); props.reset(new Opm::IncompPropertiesFromDeck(deck, global_cell)); // props.reset(new AdHocProps(param, grid->c_grid()->dimensions, grid->c_grid()->number_of_cells)); // Wells init. wells.reset(new Opm::WellsManager(deck, *grid->c_grid(), props->permeability())); // Timer init. if (deck.hasField("TSTEP")) { simtimer.init(deck); } else { simtimer.init(param); } // Gravity. gravity[2] = deck.hasField("NOGRAV") ? 0.0 : Opm::unit::gravity; // Init state variables (saturation and pressure). initStateTwophaseFromDeck(*grid->c_grid(), *props, deck, gravity[2], state); // Init polymer properties. polyprop.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 Opm::GridManager(nx, ny, nz, dx, dy, dz)); // Rock and fluid init. // props.reset(new Opm::IncompPropertiesBasic(param, grid->c_grid()->dimensions, grid->c_grid()->number_of_cells)); props.reset(new AdHocProps(param, grid->c_grid()->dimensions, grid->c_grid()->number_of_cells)); // Wells init. wells.reset(new Opm::WellsManager()); // Timer init. simtimer.init(param); // Gravity. gravity[2] = param.getDefault("gravity", 0.0); // Init state variables (saturation and pressure). initStateTwophaseBasic(*grid->c_grid(), *props, param, gravity[2], state); // 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("ads_index", Opm::PolymerProperties::NoDesorption); std::vector c_vals_visc(2, -1e100); c_vals_visc[0] = 0.0; c_vals_visc[1] = 7.0; std::vector visc_mult_vals(2, -1e100); visc_mult_vals[0] = 1.0; // polyprop.visc_mult_vals[1] = param.getDefault("c_max_viscmult", 30.0); visc_mult_vals[1] = 20.0; std::vector 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 ads_vals(3, -1e100); ads_vals[0] = 0.0; // polyprop.ads_vals[1] = param.getDefault("c_max_ads", 0.0025); ads_vals[1] = 0.0015; ads_vals[2] = 0.0025; polyprop.set(c_max, mix_param, rock_density, dead_pore_vol, res_factor, c_max_ads, static_cast(ads_index), c_vals_visc, visc_mult_vals, c_vals_ads, ads_vals); } // Initialize polymer inflow function. double poly_start = param.getDefault("poly_start_days", 300.0)*Opm::unit::day; double poly_end = param.getDefault("poly_end_days", 800.0)*Opm::unit::day; double poly_amount = param.getDefault("poly_amount", polyprop.cMax()); PolymerInflow poly_inflow(poly_start, poly_end, poly_amount); // Extra fluid init for transport solver. TwophaseFluidPolymer fluid(*props, polyprop); // 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; } } bool use_segregation_split = false; bool use_column_solver = false; bool use_gauss_seidel_gravity = false; if (use_gravity && use_reorder) { use_segregation_split = param.getDefault("use_segregation_split", use_segregation_split); if (use_segregation_split) { use_column_solver = param.getDefault("use_column_solver", use_column_solver); if (use_column_solver) { // use_gauss_seidel_gravity is not implemented for polymer. use_gauss_seidel_gravity = param.getDefault("use_gauss_seidel_gravity", use_gauss_seidel_gravity); } } } // Source-related variables init. int num_cells = grid->c_grid()->number_of_cells; std::vector totmob; std::vector omega; // Will remain empty if no gravity. // Extra rock init. std::vector porevol; computePorevolume(*grid->c_grid(), *props, porevol); double tot_porevol_init = std::accumulate(porevol.begin(), porevol.end(), 0.0); // We need a separate reorder_sat, because the reorder // code expects a scalar sw, not both sw and so. std::vector reorder_sat(num_cells); std::vector src(num_cells, 0.0); // Initialising src if (wells->c_wells()) { Opm::wellsToSrc(*wells->c_wells(), num_cells, src); } else { const double default_injection = use_gravity ? 0.0 : 0.1; const double flow_per_sec = param.getDefault("injected_porevolumes_per_day", default_injection) *tot_porevol_init/Opm::unit::day; src[0] = flow_per_sec; src[num_cells - 1] = -flow_per_sec; } TransportSource* tsrc = create_transport_source(2, 2); double ssrc[] = { 1.0, 0.0 }; double ssink[] = { 0.0, 1.0 }; double zdummy[] = { 0.0, 0.0 }; for (int cell = 0; cell < num_cells; ++cell) { if (src[cell] > 0.0) { append_transport_source(cell, 2, 0, src[cell], ssrc, zdummy, tsrc); } else if (src[cell] < 0.0) { append_transport_source(cell, 2, 0, src[cell], ssink, zdummy, tsrc); } } std::vector reorder_src = src; // Boundary conditions. Opm::FlowBCManager bcs; if (param.getDefault("use_pside", false)) { int pside = param.get("pside"); double pside_pressure = param.get("pside_pressure"); bcs.pressureSide(*grid->c_grid(), Opm::FlowBCManager::Side(pside), pside_pressure); } // Solvers init. // Pressure solver. #ifdef EXPERIMENT_ISTL Opm::LinearSolverIstl linsolver(param); #else Opm::LinearSolverUmfpack linsolver; #endif // EXPERIMENT_ISTL const double *grav = use_gravity ? &gravity[0] : 0; Opm::IncompTpfa psolver(*grid->c_grid(), props->permeability(), grav, linsolver); // Reordering solver. const double nl_tolerance = param.getDefault("nl_tolerance", 1e-9); const int nl_maxiter = param.getDefault("nl_maxiter", 30); Opm::TransportModelPolymer::SingleCellMethod method; std::string method_string = param.getDefault("single_cell_method", std::string("Bracketing")); if (method_string == "Bracketing") { method = Opm::TransportModelPolymer::Bracketing; } else if (method_string == "Newton") { method = Opm::TransportModelPolymer::Newton; } else { THROW("Unknown method: " << method_string); } Opm::TransportModelPolymer reorder_model(*grid->c_grid(), props->porosity(), &porevol[0], *props, polyprop, method, nl_tolerance, nl_maxiter); if (use_gauss_seidel_gravity) { THROW("use_gauss_seidel_gravity option not yet implemented for polymer case."); // reorder_model.initGravity(grav); } // Non-reordering solver. NewtonPolymerTransportModel model (fluid, *grid->c_grid(), porevol, grav, guess_old_solution); if (use_gravity) { model.initGravityTrans(*grid->c_grid(), psolver.getHalfTrans()); } TransportSolver tsolver(model); // Column-based gravity segregation solver. typedef std::pair, std::vector > > ColMap; ColMap columns; if (use_column_solver) { Opm::extractColumn(*grid->c_grid(), columns); } Opm::GravityColumnSolverPolymer colsolver(model, *grid->c_grid(), nl_tolerance, nl_maxiter); // // // Not implemented for polymer. // // Control init. // Opm::ImplicitTransportDetails::NRReport rpt; // Opm::ImplicitTransportDetails::NRControl ctrl; // if (!use_reorder || use_segregation_split) { // ctrl.max_it = param.getDefault("max_it", 20); // ctrl.verbosity = param.getDefault("verbosity", 0); // ctrl.max_it_ls = param.getDefault("max_it_ls", 5); // } // // Linear solver init. // using Opm::ImplicitTransportLinAlgSupport::CSRMatrixUmfpackSolver; // CSRMatrixUmfpackSolver linsolve; // The allcells vector is used in calls to computeTotalMobility() // and computeTotalMobilityOmega(). std::vector allcells(num_cells); for (int cell = 0; cell < num_cells; ++cell) { allcells[cell] = cell; } // Warn if any parameters are unused. if (param.anyUnused()) { std::cout << "-------------------- Unused parameters: --------------------\n"; param.displayUsage(); std::cout << "----------------------------------------------------------------" << std::endl; } // Write parameters used for later reference. if (output) { param.writeParam(output_dir + "/spu_2p.param"); } // 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(); std::cout << "\n\n================ Starting main simulation loop ===============" << std::endl; 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); for (; !simtimer.done(); ++simtimer) { // Report timestep and (optionally) write state to disk. simtimer.report(std::cout); if (output && (simtimer.currentStepNum() % output_interval == 0)) { outputState(*grid->c_grid(), state, simtimer.currentStepNum(), output_dir); } // Solve pressure. if (use_gravity) { computeTotalMobilityOmega(*props, polyprop, allcells, state.saturation(), state.concentration(), totmob, omega); } else { computeTotalMobility(*props, polyprop, allcells, state.saturation(), state.concentration(), totmob); } pressure_timer.start(); psolver.solve(totmob, omega, src, bcs.c_bcs(), state.pressure(), state.faceflux()); pressure_timer.stop(); double pt = pressure_timer.secsSinceStart(); std::cout << "Pressure solver took: " << pt << " seconds." << std::endl; ptime += pt; // Process transport sources (to include bdy terms). if (use_reorder) { Opm::computeTransportSource(*grid->c_grid(), src, state.faceflux(), 1.0, reorder_src); } else { clear_transport_source(tsrc); for (int cell = 0; cell < num_cells; ++cell) { if (src[cell] > 0.0) { append_transport_source(cell, 2, 0, src[cell], ssrc, zdummy, tsrc); } else if (src[cell] < 0.0) { append_transport_source(cell, 2, 0, src[cell], ssink, zdummy, tsrc); } } } // Find inflow rate. const double current_time = simtimer.currentTime(); const double stepsize = simtimer.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 (use_reorder) { Opm::toWaterSat(state.saturation(), reorder_sat); reorder_model.solve(&state.faceflux()[0], &reorder_src[0], stepsize, inflow_c, &reorder_sat[0], &state.concentration()[0], &state.maxconcentration()[0]); Opm::toBothSat(reorder_sat, state.saturation()); Opm::computeInjectedProduced(*props, state.saturation(), src, simtimer.currentStepLength(), injected, produced); if (use_segregation_split) { if (use_column_solver) { if (use_gauss_seidel_gravity) { THROW("use_gauss_seidel_gravity option not implemented for polymer."); // reorder_model.solveGravity(columns, simtimer.currentStepLength(), reorder_sat); // Opm::toBothSat(reorder_sat, state.saturation()); } else { colsolver.solve(columns, simtimer.currentStepLength(), state.saturation(), state.concentration(), state.maxconcentration()); } } else { THROW("use_segregation_split option for polymer is only implemented in the use_column_solver case."); // std::vector fluxes = state.faceflux(); // std::fill(state.faceflux().begin(), state.faceflux().end(), 0.0); // tsolver.solve(*grid->c_grid(), tsrc, simtimer.currentStepLength(), ctrl, state, linsolve, rpt); // std::cout << rpt; // state.faceflux() = fluxes; } } } else { THROW("Non-reordering transport solver not implemented for polymer."); // tsolver.solve(*grid->c_grid(), tsrc, simtimer.currentStepLength(), ctrl, state, linsolve, rpt); // std::cout << rpt; // Opm::computeInjectedProduced(*props, state.saturation(), src, simtimer.currentStepLength(), injected, produced); } 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(), polyprop.deadPoreVol()); polymass_adsorbed = Opm::computePolymerAdsorbed(*props, polyprop, porevol, state.maxconcentration()); Opm::computeInjectedProduced(*props, polyprop, state.saturation(), state.concentration(), src, simtimer.currentStepLength(), inflow_c, injected, produced, polyinj, polyprod); 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(simtimer.currentTime() + simtimer.currentStepLength(), produced[0]/(produced[0] + produced[1]), tot_produced[0]/tot_porevol_init); } total_timer.stop(); std::cout << "\n\n================ End of simulation ===============\n" << "Total time taken: " << total_timer.secsSinceStart() << "\n Pressure time: " << ptime << "\n Transport time: " << ttime << std::endl; if (output) { outputState(*grid->c_grid(), state, simtimer.currentStepNum(), output_dir); outputWaterCut(watercut, output_dir); } destroy_transport_source(tsrc); }