2011-10-05 20:06:03 +02:00
|
|
|
/*===========================================================================
|
|
|
|
|
//
|
|
|
|
|
// File: spu_2p.cpp
|
|
|
|
|
//
|
|
|
|
|
// Created: 2011-10-05 10:29:01+0200
|
|
|
|
|
//
|
|
|
|
|
// Authors: Ingeborg S. Ligaarden <Ingeborg.Ligaarden@sintef.no>
|
|
|
|
|
// Jostein R. Natvig <Jostein.R.Natvig@sintef.no>
|
|
|
|
|
// Halvor M. Nilsen <HalvorMoll.Nilsen@sintef.no>
|
|
|
|
|
// Atgeirr F. Rasmussen <atgeirr@sintef.no>
|
|
|
|
|
// Bård Skaflestad <Bard.Skaflestad@sintef.no>
|
|
|
|
|
//
|
|
|
|
|
//==========================================================================*/
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/*
|
|
|
|
|
Copyright 2011 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
Copyright 2011 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 <http://www.gnu.org/licenses/>.
|
2012-03-06 16:37:49 +01:00
|
|
|
*/
|
2012-02-21 14:46:28 +01:00
|
|
|
|
|
|
|
|
#if HAVE_CONFIG_H
|
2011-12-05 19:22:20 +01:00
|
|
|
#include "config.h"
|
2012-02-21 14:46:28 +01:00
|
|
|
#endif // HAVE_CONFIG_H
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2011-12-12 11:28:09 +01:00
|
|
|
#include <opm/core/linalg/sparse_sys.h>
|
2012-02-21 21:45:04 +01:00
|
|
|
#include <opm/core/linalg/LinearSolverUmfpack.hpp>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-02-21 22:14:01 +01:00
|
|
|
// #define EXPERIMENT_ISTL
|
|
|
|
|
#ifdef EXPERIMENT_ISTL
|
|
|
|
|
#include <opm/core/linalg/LinearSolverIstl.hpp>
|
|
|
|
|
#endif
|
|
|
|
|
|
2012-02-20 13:42:42 +01:00
|
|
|
#include <opm/core/pressure/IncompTpfa.hpp>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-02-20 11:34:51 +01:00
|
|
|
#include <opm/core/GridManager.hpp>
|
|
|
|
|
#include <opm/core/grid.h>
|
2012-01-19 16:48:27 +01:00
|
|
|
#include <opm/core/utility/ErrorMacros.hpp>
|
2012-02-09 23:15:14 +01:00
|
|
|
#include <opm/core/utility/StopWatch.hpp>
|
2012-01-19 14:12:53 +01:00
|
|
|
#include <opm/core/utility/Units.hpp>
|
2012-02-19 21:32:35 +01:00
|
|
|
#include <opm/core/utility/writeVtkData.hpp>
|
2012-02-26 21:05:19 +01:00
|
|
|
#include <opm/core/utility/miscUtilities.hpp>
|
2012-01-19 13:58:25 +01:00
|
|
|
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2011-12-12 11:28:09 +01:00
|
|
|
#include <opm/core/fluid/SimpleFluid2p.hpp>
|
2012-01-23 13:47:00 +01:00
|
|
|
#include <opm/core/fluid/IncompPropertiesBasic.hpp>
|
2012-01-31 09:44:52 +01:00
|
|
|
#include <opm/core/fluid/IncompPropertiesFromDeck.hpp>
|
2012-03-06 16:37:49 +01:00
|
|
|
|
2011-12-12 11:28:09 +01:00
|
|
|
#include <opm/core/transport/transport_source.h>
|
|
|
|
|
#include <opm/core/transport/CSRMatrixUmfpackSolver.hpp>
|
|
|
|
|
#include <opm/core/transport/NormSupport.hpp>
|
|
|
|
|
#include <opm/core/transport/ImplicitAssembly.hpp>
|
|
|
|
|
#include <opm/core/transport/ImplicitTransport.hpp>
|
|
|
|
|
#include <opm/core/transport/JacobianSystem.hpp>
|
|
|
|
|
#include <opm/core/transport/CSRMatrixBlockAssembler.hpp>
|
|
|
|
|
#include <opm/core/transport/SinglePointUpwindTwoPhase.hpp>
|
2011-10-13 19:38:33 +02:00
|
|
|
|
2012-03-02 13:55:54 +01:00
|
|
|
#include <opm/core/ColumnExtract.hpp>
|
|
|
|
|
#include <opm/core/transport/GravityColumnSolver.hpp>
|
|
|
|
|
|
2012-02-10 10:48:18 +01:00
|
|
|
#include <opm/core/transport/reorder/TransportModelTwophase.hpp>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-24 15:22:51 +01:00
|
|
|
#include <boost/filesystem/convenience.hpp>
|
2012-01-27 11:43:14 +01:00
|
|
|
#include <boost/scoped_ptr.hpp>
|
2012-01-31 22:33:48 +01:00
|
|
|
#include <boost/lexical_cast.hpp>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-24 15:22:51 +01:00
|
|
|
#include <cassert>
|
|
|
|
|
#include <cstddef>
|
|
|
|
|
|
|
|
|
|
#include <algorithm>
|
|
|
|
|
#include <tr1/array>
|
|
|
|
|
#include <functional>
|
|
|
|
|
#include <iostream>
|
|
|
|
|
#include <iomanip>
|
|
|
|
|
#include <fstream>
|
|
|
|
|
#include <iterator>
|
|
|
|
|
#include <vector>
|
2011-10-05 20:06:03 +02:00
|
|
|
|
|
|
|
|
|
2012-01-31 09:44:52 +01:00
|
|
|
|
|
|
|
|
|
2011-10-05 20:06:03 +02:00
|
|
|
class ReservoirState {
|
|
|
|
|
public:
|
2012-01-19 16:48:27 +01:00
|
|
|
ReservoirState(const UnstructuredGrid* g, const int num_phases = 2)
|
2012-03-06 16:37:49 +01:00
|
|
|
: press_ (g->number_of_cells, 0.0),
|
|
|
|
|
fpress_(g->number_of_faces, 0.0),
|
|
|
|
|
flux_ (g->number_of_faces, 0.0),
|
|
|
|
|
sat_ (num_phases * g->number_of_cells, 0.0)
|
2012-02-02 14:24:49 +01:00
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
for (int cell = 0; cell < g->number_of_cells; ++cell) {
|
|
|
|
|
sat_[num_phases*cell + num_phases - 1] = 1.0;
|
|
|
|
|
}
|
2012-02-02 14:24:49 +01:00
|
|
|
}
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-19 16:48:27 +01:00
|
|
|
int numPhases() const { return sat_.size()/press_.size(); }
|
|
|
|
|
|
2011-10-05 20:06:03 +02:00
|
|
|
::std::vector<double>& pressure () { return press_ ; }
|
|
|
|
|
::std::vector<double>& facepressure() { return fpress_; }
|
|
|
|
|
::std::vector<double>& faceflux () { return flux_ ; }
|
|
|
|
|
::std::vector<double>& saturation () { return sat_ ; }
|
|
|
|
|
|
2012-01-19 16:48:27 +01:00
|
|
|
const ::std::vector<double>& pressure () const { return press_ ; }
|
|
|
|
|
const ::std::vector<double>& facepressure() const { return fpress_; }
|
|
|
|
|
const ::std::vector<double>& faceflux () const { return flux_ ; }
|
|
|
|
|
const ::std::vector<double>& saturation () const { return sat_ ; }
|
2011-10-05 20:06:03 +02:00
|
|
|
|
|
|
|
|
private:
|
|
|
|
|
::std::vector<double> press_ ;
|
|
|
|
|
::std::vector<double> fpress_;
|
|
|
|
|
::std::vector<double> flux_ ;
|
|
|
|
|
::std::vector<double> sat_ ;
|
|
|
|
|
};
|
|
|
|
|
|
2012-01-23 13:47:00 +01:00
|
|
|
|
|
|
|
|
|
2012-01-19 13:58:25 +01:00
|
|
|
template <class State>
|
2012-01-31 22:33:48 +01:00
|
|
|
void outputState(const UnstructuredGrid* grid,
|
2012-03-06 16:37:49 +01:00
|
|
|
const State& state,
|
|
|
|
|
const int step,
|
|
|
|
|
const std::string& output_dir)
|
2012-01-19 16:48:27 +01:00
|
|
|
{
|
2012-02-19 21:32:35 +01:00
|
|
|
// Write data in VTK format.
|
2012-01-19 16:48:27 +01:00
|
|
|
std::ostringstream vtkfilename;
|
2012-01-31 22:33:48 +01:00
|
|
|
vtkfilename << output_dir << "/output-" << std::setw(3) << std::setfill('0') << step << ".vtu";
|
2012-01-19 16:48:27 +01:00
|
|
|
std::ofstream vtkfile(vtkfilename.str().c_str());
|
|
|
|
|
if (!vtkfile) {
|
2012-03-06 16:37:49 +01:00
|
|
|
THROW("Failed to open " << vtkfilename.str());
|
2012-01-19 16:48:27 +01:00
|
|
|
}
|
2012-02-19 21:32:35 +01:00
|
|
|
Opm::DataMap dm;
|
|
|
|
|
dm["saturation"] = &state.saturation();
|
|
|
|
|
dm["pressure"] = &state.pressure();
|
2012-02-25 22:29:35 +01:00
|
|
|
std::vector<double> cell_velocity;
|
2012-02-26 21:05:19 +01:00
|
|
|
Opm::estimateCellVelocity(*grid, state.faceflux(), cell_velocity);
|
2012-02-25 22:29:35 +01:00
|
|
|
dm["velocity"] = &cell_velocity;
|
2012-02-19 21:32:35 +01:00
|
|
|
Opm::writeVtkData(grid, dm, vtkfile);
|
|
|
|
|
|
|
|
|
|
// Write data (not grid) in Matlab format
|
|
|
|
|
for (Opm::DataMap::const_iterator it = dm.begin(); it != dm.end(); ++it) {
|
2012-03-06 16:37:49 +01:00
|
|
|
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<double>& d = *(it->second);
|
|
|
|
|
std::copy(d.begin(), d.end(), std::ostream_iterator<double>(file, "\n"));
|
2012-01-31 22:33:48 +01:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2012-01-19 13:58:25 +01:00
|
|
|
|
2012-01-27 09:19:59 +01:00
|
|
|
// --------------- Types needed to define transport solver ---------------
|
|
|
|
|
|
2012-01-27 11:43:14 +01:00
|
|
|
class SimpleFluid2pWrappingProps
|
|
|
|
|
{
|
|
|
|
|
public:
|
|
|
|
|
SimpleFluid2pWrappingProps(const Opm::IncompPropertiesInterface& props)
|
2012-03-06 16:37:49 +01:00
|
|
|
: props_(props),
|
|
|
|
|
smin_(props.numCells()*props.numPhases()),
|
|
|
|
|
smax_(props.numCells()*props.numPhases())
|
2012-01-27 11:43:14 +01:00
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
if (props.numPhases() != 2) {
|
|
|
|
|
THROW("SimpleFluid2pWrapper requires 2 phases.");
|
|
|
|
|
}
|
|
|
|
|
const int num_cells = props.numCells();
|
|
|
|
|
std::vector<int> cells(num_cells);
|
|
|
|
|
for (int c = 0; c < num_cells; ++c) {
|
|
|
|
|
cells[c] = c;
|
|
|
|
|
}
|
|
|
|
|
props.satRange(num_cells, &cells[0], &smin_[0], &smax_[0]);
|
2012-01-27 11:43:14 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
double density(int phase) const
|
|
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
return props_.density()[phase];
|
2012-01-27 11:43:14 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
template <class Sat,
|
2012-03-06 16:37:49 +01:00
|
|
|
class Mob,
|
|
|
|
|
class DMob>
|
2012-01-27 11:43:14 +01:00
|
|
|
void mobility(int c, const Sat& s, Mob& mob, DMob& dmob) const
|
|
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
props_.relperm(1, &s[0], &c, &mob[0], &dmob[0]);
|
|
|
|
|
const double* mu = props_.viscosity();
|
|
|
|
|
mob[0] /= mu[0];
|
|
|
|
|
mob[1] /= mu[1];
|
|
|
|
|
// Recall that we use Fortran ordering for kr derivatives,
|
|
|
|
|
// therefore dmob[i*2 + j] is row j and column i of the
|
|
|
|
|
// matrix.
|
|
|
|
|
// Each row corresponds to a kr function, so which mu to
|
|
|
|
|
// divide by also depends on the row, j.
|
|
|
|
|
dmob[0*2 + 0] /= mu[0];
|
|
|
|
|
dmob[0*2 + 1] /= mu[1];
|
|
|
|
|
dmob[1*2 + 0] /= mu[0];
|
|
|
|
|
dmob[1*2 + 1] /= mu[1];
|
2012-01-27 11:43:14 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
template <class Sat,
|
2012-03-06 16:37:49 +01:00
|
|
|
class Pcap,
|
|
|
|
|
class DPcap>
|
2012-01-27 11:43:14 +01:00
|
|
|
void pc(int c, const Sat& s, Pcap& pcap, DPcap& dpcap) const
|
|
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
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);
|
2012-01-27 11:43:14 +01:00
|
|
|
}
|
2012-03-06 16:37:49 +01:00
|
|
|
|
2012-02-20 13:50:45 +01:00
|
|
|
double s_min(int c) const
|
|
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
return smin_[c];
|
2012-02-20 13:50:45 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
double s_max(int c) const
|
|
|
|
|
{
|
2012-03-06 16:37:49 +01:00
|
|
|
return smax_[c];
|
2012-02-20 13:50:45 +01:00
|
|
|
}
|
2012-01-27 11:43:14 +01:00
|
|
|
|
|
|
|
|
private:
|
|
|
|
|
const Opm::IncompPropertiesInterface& props_;
|
2012-02-20 13:50:45 +01:00
|
|
|
std::vector<double> smin_;
|
|
|
|
|
std::vector<double> smax_;
|
2012-01-27 11:43:14 +01:00
|
|
|
};
|
|
|
|
|
|
|
|
|
|
typedef SimpleFluid2pWrappingProps TwophaseFluid;
|
2012-01-27 09:19:59 +01:00
|
|
|
typedef Opm::SinglePointUpwindTwoPhase<TwophaseFluid> TransportModel;
|
|
|
|
|
|
|
|
|
|
using namespace Opm::ImplicitTransportDefault;
|
|
|
|
|
|
|
|
|
|
typedef NewtonVectorCollection< ::std::vector<double> > NVecColl;
|
|
|
|
|
typedef JacobianSystem < struct CSRMatrix, NVecColl > JacSys;
|
|
|
|
|
|
|
|
|
|
template <class Vector>
|
|
|
|
|
class MaxNorm {
|
|
|
|
|
public:
|
|
|
|
|
static double
|
|
|
|
|
norm(const Vector& v) {
|
|
|
|
|
return AccumulationNorm <Vector, MaxAbs>::norm(v);
|
|
|
|
|
}
|
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
typedef Opm::ImplicitTransport<TransportModel,
|
2012-03-06 16:37:49 +01:00
|
|
|
JacSys ,
|
|
|
|
|
MaxNorm ,
|
|
|
|
|
VectorNegater ,
|
|
|
|
|
VectorZero ,
|
|
|
|
|
MatrixZero ,
|
|
|
|
|
VectorAssign > TransportSolver;
|
2012-01-27 09:19:59 +01:00
|
|
|
|
|
|
|
|
|
2012-01-23 13:47:00 +01:00
|
|
|
|
2012-01-19 13:58:25 +01:00
|
|
|
// ----------------- Main program -----------------
|
2011-10-05 20:06:03 +02:00
|
|
|
int
|
2012-01-19 13:58:25 +01:00
|
|
|
main(int argc, char** argv)
|
2011-10-05 20:06:03 +02:00
|
|
|
{
|
2012-01-27 09:19:59 +01:00
|
|
|
std::cout << "\n================ Test program for incompressible two-phase flow ===============\n\n";
|
2012-01-19 14:12:53 +01:00
|
|
|
Opm::parameter::ParameterGroup param(argc, argv, false);
|
2012-01-27 09:19:59 +01:00
|
|
|
std::cout << "--------------- Reading parameters ---------------" << std::endl;
|
2012-01-27 11:43:14 +01:00
|
|
|
|
|
|
|
|
// Reading various control parameters.
|
2012-01-19 13:58:25 +01:00
|
|
|
const int num_psteps = param.getDefault("num_psteps", 1);
|
2012-01-24 10:07:10 +01:00
|
|
|
const double stepsize_days = param.getDefault("stepsize_days", 1.0);
|
2012-01-19 23:51:09 +01:00
|
|
|
const double stepsize = Opm::unit::convert::from(stepsize_days, Opm::unit::day);
|
|
|
|
|
const bool guess_old_solution = param.getDefault("guess_old_solution", false);
|
2012-01-31 22:35:50 +01:00
|
|
|
const bool use_reorder = param.getDefault("use_reorder", true);
|
2012-01-20 15:51:29 +01:00
|
|
|
const bool output = param.getDefault("output", true);
|
2012-01-24 15:22:51 +01:00
|
|
|
std::string output_dir;
|
|
|
|
|
if (output) {
|
2012-03-06 16:37:49 +01:00
|
|
|
output_dir = param.getDefault("output_dir", std::string("output"));
|
|
|
|
|
// Ensure that output dir exists
|
|
|
|
|
boost::filesystem::path fpath(output_dir);
|
|
|
|
|
create_directories(fpath);
|
2012-01-24 15:22:51 +01:00
|
|
|
}
|
2012-01-19 16:48:27 +01:00
|
|
|
|
2012-01-27 11:43:14 +01:00
|
|
|
// If we have a "deck_filename", grid and props will be read from that.
|
2012-01-31 09:44:52 +01:00
|
|
|
bool use_deck = param.has("deck_filename");
|
2012-02-20 11:34:51 +01:00
|
|
|
boost::scoped_ptr<Opm::GridManager> grid;
|
2012-01-27 11:43:14 +01:00
|
|
|
boost::scoped_ptr<Opm::IncompPropertiesInterface> props;
|
2012-01-31 09:44:52 +01:00
|
|
|
if (use_deck) {
|
2012-03-06 16:37:49 +01:00
|
|
|
std::string deck_filename = param.get<std::string>("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<int> global_cell(gc, gc + grid->c_grid()->number_of_cells);
|
|
|
|
|
props.reset(new Opm::IncompPropertiesFromDeck(deck, global_cell));
|
2012-01-31 09:44:52 +01:00
|
|
|
} else {
|
2012-03-06 16:37:49 +01:00
|
|
|
// 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));
|
2012-01-31 09:44:52 +01:00
|
|
|
}
|
2012-01-27 11:43:14 +01:00
|
|
|
|
|
|
|
|
// Extra rock init.
|
2012-01-23 13:06:25 +01:00
|
|
|
std::vector<double> porevol;
|
2012-02-26 21:05:19 +01:00
|
|
|
computePorevolume(*grid->c_grid(), *props, porevol);
|
2012-01-24 09:52:15 +01:00
|
|
|
double tot_porevol = std::accumulate(porevol.begin(), porevol.end(), 0.0);
|
2012-01-23 13:06:25 +01:00
|
|
|
|
2012-01-31 09:44:52 +01:00
|
|
|
// Extra fluid init for transport solver.
|
2012-01-27 11:43:14 +01:00
|
|
|
TwophaseFluid fluid(*props);
|
2011-10-06 17:25:13 +02:00
|
|
|
|
2012-02-15 22:43:56 +01:00
|
|
|
// Gravity init.
|
|
|
|
|
double gravity[3] = { 0.0 };
|
|
|
|
|
double g = param.getDefault("gravity", 0.0);
|
|
|
|
|
bool use_gravity = g != 0.0;
|
|
|
|
|
if (use_gravity) {
|
2012-03-06 16:37:49 +01:00
|
|
|
gravity[grid->c_grid()->dimensions - 1] = g;
|
|
|
|
|
if (props->density()[0] == props->density()[1]) {
|
|
|
|
|
std::cout << "**** Warning: nonzero gravity, but zero density difference." << std::endl;
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
}
|
2012-03-01 10:32:58 +01:00
|
|
|
bool use_segregation_split = false;
|
2012-03-02 13:55:54 +01:00
|
|
|
bool use_column_solver = false;
|
2012-03-01 10:32:58 +01:00
|
|
|
if (use_gravity && use_reorder) {
|
2012-03-06 16:37:49 +01:00
|
|
|
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);
|
|
|
|
|
}
|
2012-03-01 10:32:58 +01:00
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
|
2012-01-23 13:06:25 +01:00
|
|
|
// Solvers init.
|
2012-02-10 10:48:18 +01:00
|
|
|
// Pressure solver.
|
2012-02-21 22:14:01 +01:00
|
|
|
#ifdef EXPERIMENT_ISTL
|
|
|
|
|
Opm::LinearSolverIstl linsolver(param);
|
|
|
|
|
#else
|
2012-02-21 21:45:04 +01:00
|
|
|
Opm::LinearSolverUmfpack linsolver;
|
2012-02-21 22:14:01 +01:00
|
|
|
#endif // EXPERIMENT_ISTL
|
2012-02-27 17:34:57 +01:00
|
|
|
const double *grav = use_gravity ? &gravity[0] : 0;
|
|
|
|
|
Opm::IncompTpfa psolver(*grid->c_grid(), props->permeability(), grav, linsolver);
|
2012-02-10 10:48:18 +01:00
|
|
|
// Non-reordering solver.
|
2012-02-27 17:34:57 +01:00
|
|
|
TransportModel model (fluid, *grid->c_grid(), porevol, grav, guess_old_solution);
|
|
|
|
|
if (use_gravity) {
|
|
|
|
|
model.initGravityTrans(*grid->c_grid(), psolver.getHalfTrans());
|
|
|
|
|
}
|
2012-01-23 13:06:25 +01:00
|
|
|
TransportSolver tsolver(model);
|
2012-02-10 10:48:18 +01:00
|
|
|
// Reordering solver.
|
2012-02-17 14:18:03 +01:00
|
|
|
const double nltol = param.getDefault("nl_tolerance", 1e-9);
|
|
|
|
|
const int maxit = param.getDefault("nl_maxiter", 30);
|
|
|
|
|
Opm::TransportModelTwophase reorder_model(*grid->c_grid(), &porevol[0], *props, nltol, maxit);
|
2012-03-02 13:55:54 +01:00
|
|
|
// Column-based gravity segregation solver.
|
|
|
|
|
typedef std::map<int, std::vector<int> > ColMap;
|
|
|
|
|
ColMap columns;
|
|
|
|
|
if (use_column_solver) {
|
2012-03-06 16:37:49 +01:00
|
|
|
Opm::extractColumn(*grid->c_grid(), columns);
|
2012-03-02 13:55:54 +01:00
|
|
|
}
|
2012-03-06 22:33:19 +01:00
|
|
|
Opm::GravityColumnSolver<TransportModel> colsolver(model, *grid->c_grid(), nltol, maxit);
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-23 13:06:25 +01:00
|
|
|
// State-related and source-related variables init.
|
2012-02-15 22:43:56 +01:00
|
|
|
int num_cells = grid->c_grid()->number_of_cells;
|
2012-02-02 14:24:49 +01:00
|
|
|
std::vector<double> totmob;
|
2012-02-20 13:42:42 +01:00
|
|
|
std::vector<double> omega; // Will remain empty if no gravity.
|
2012-02-02 14:24:49 +01:00
|
|
|
ReservoirState state(grid->c_grid(), props->numPhases());
|
2012-01-23 13:06:25 +01:00
|
|
|
// We need a separate reorder_sat, because the reorder
|
2012-01-23 12:43:37 +01:00
|
|
|
// code expects a scalar sw, not both sw and so.
|
2012-02-15 22:43:56 +01:00
|
|
|
std::vector<double> reorder_sat(num_cells);
|
|
|
|
|
std::vector<double> src(num_cells, 0.0);
|
|
|
|
|
int scenario = param.getDefault("scenario", 0);
|
|
|
|
|
switch (scenario) {
|
|
|
|
|
case 0:
|
2012-03-06 22:33:19 +01:00
|
|
|
{
|
|
|
|
|
std::cout << "==== Scenario 0: single-cell source and sink.\n";
|
|
|
|
|
double flow_per_sec = 0.1*tot_porevol/Opm::unit::day;
|
|
|
|
|
src[0] = flow_per_sec;
|
|
|
|
|
src[grid->c_grid()->number_of_cells - 1] = -flow_per_sec;
|
|
|
|
|
break;
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
case 1:
|
2012-03-06 22:33:19 +01:00
|
|
|
{
|
|
|
|
|
std::cout << "==== Scenario 1: half source, half sink.\n";
|
|
|
|
|
double flow_per_sec = 0.1*porevol[0]/Opm::unit::day;
|
|
|
|
|
std::fill(src.begin(), src.begin() + src.size()/2, flow_per_sec);
|
|
|
|
|
std::fill(src.begin() + src.size()/2, src.end(), -flow_per_sec);
|
|
|
|
|
break;
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
case 2:
|
2012-03-06 22:33:19 +01:00
|
|
|
{
|
|
|
|
|
std::cout << "==== Scenario 2: gravity convection.\n";
|
|
|
|
|
if (!use_gravity) {
|
|
|
|
|
std::cout << "**** Warning: running gravity convection scenario, but gravity is zero." << std::endl;
|
|
|
|
|
}
|
|
|
|
|
if (use_deck) {
|
|
|
|
|
std::cout << "**** Warning: running gravity convection scenario, which expects a cartesian grid."
|
|
|
|
|
<< std::endl;
|
|
|
|
|
}
|
|
|
|
|
std::vector<double>& sat = state.saturation();
|
|
|
|
|
const int *glob_cell = grid->c_grid()->global_cell;
|
|
|
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
|
|
|
|
const int* cd = grid->c_grid()->cartdims;
|
|
|
|
|
const int gc = glob_cell == 0 ? cell : glob_cell[cell];
|
|
|
|
|
bool left = (gc % cd[0]) < cd[0]/2;
|
|
|
|
|
sat[2*cell] = left ? 1.0 : 0.0;
|
|
|
|
|
sat[2*cell + 1] = 1.0 - sat[2*cell];
|
|
|
|
|
}
|
|
|
|
|
break;
|
|
|
|
|
}
|
2012-03-06 16:37:49 +01:00
|
|
|
case 3:
|
2012-03-06 22:33:19 +01:00
|
|
|
{
|
|
|
|
|
std::cout << "==== Scenario 3: gravity segregation.\n";
|
|
|
|
|
if (!use_gravity) {
|
|
|
|
|
std::cout << "**** Warning: running gravity segregation scenario, but gravity is zero." << std::endl;
|
|
|
|
|
}
|
|
|
|
|
if (use_deck) {
|
|
|
|
|
std::cout << "**** Warning: running gravity segregation scenario, which expects a cartesian grid."
|
|
|
|
|
<< std::endl;
|
|
|
|
|
}
|
|
|
|
|
std::vector<double>& sat = state.saturation();
|
|
|
|
|
const int *glob_cell = grid->c_grid()->global_cell;
|
|
|
|
|
// Water on top
|
|
|
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
|
|
|
|
const int* cd = grid->c_grid()->cartdims;
|
|
|
|
|
const int gc = glob_cell == 0 ? cell : glob_cell[cell];
|
|
|
|
|
bool top = (gc / cd[0] / cd[1]) < cd[2]/2;
|
|
|
|
|
sat[2*cell] = top ? 1.0 : 0.0;
|
|
|
|
|
sat[2*cell + 1 ] = 1.0 - sat[2*cell];
|
|
|
|
|
}
|
|
|
|
|
break;
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
default:
|
2012-03-06 22:33:19 +01:00
|
|
|
{
|
|
|
|
|
THROW("==== Scenario " << scenario << " is unknown.");
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
}
|
2011-10-05 20:06:03 +02:00
|
|
|
TransportSource* tsrc = create_transport_source(2, 2);
|
2012-01-24 09:58:24 +01:00
|
|
|
double ssrc[] = { 1.0, 0.0 };
|
|
|
|
|
double ssink[] = { 0.0, 1.0 };
|
2011-10-05 20:06:03 +02:00
|
|
|
double zdummy[] = { 0.0, 0.0 };
|
2012-02-15 22:43:56 +01:00
|
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
2012-03-06 16:37:49 +01:00
|
|
|
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);
|
|
|
|
|
}
|
2012-02-15 22:43:56 +01:00
|
|
|
}
|
2012-01-23 12:43:37 +01:00
|
|
|
std::vector<double> reorder_src = src;
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-23 13:06:25 +01:00
|
|
|
// Control init.
|
2012-01-19 13:58:25 +01:00
|
|
|
Opm::ImplicitTransportDetails::NRReport rpt;
|
|
|
|
|
Opm::ImplicitTransportDetails::NRControl ctrl;
|
2012-01-19 14:12:53 +01:00
|
|
|
double current_time = 0.0;
|
|
|
|
|
double total_time = stepsize*num_psteps;
|
2012-03-01 10:32:58 +01:00
|
|
|
if (!use_reorder || use_segregation_split) {
|
2012-03-06 16:37:49 +01:00
|
|
|
ctrl.max_it = param.getDefault("max_it", 20);
|
|
|
|
|
ctrl.verbosity = param.getDefault("verbosity", 0);
|
|
|
|
|
ctrl.max_it_ls = param.getDefault("max_it_ls", 5);
|
2012-02-17 14:18:03 +01:00
|
|
|
}
|
2011-10-05 20:06:03 +02:00
|
|
|
|
2012-01-23 13:06:25 +01:00
|
|
|
// Linear solver init.
|
2012-01-19 13:58:25 +01:00
|
|
|
using Opm::ImplicitTransportLinAlgSupport::CSRMatrixUmfpackSolver;
|
|
|
|
|
CSRMatrixUmfpackSolver linsolve;
|
|
|
|
|
|
2012-02-26 21:05:19 +01:00
|
|
|
// The allcells vector is used in calls to computeTotalMobility()
|
|
|
|
|
// and computeTotalMobilityOmega().
|
|
|
|
|
std::vector<int> allcells(num_cells);
|
|
|
|
|
for (int cell = 0; cell < num_cells; ++cell) {
|
2012-03-06 16:37:49 +01:00
|
|
|
allcells[cell] = cell;
|
2012-02-26 21:05:19 +01:00
|
|
|
}
|
|
|
|
|
|
2012-01-24 16:00:16 +01:00
|
|
|
// Warn if any parameters are unused.
|
|
|
|
|
if (param.anyUnused()) {
|
2012-03-06 16:37:49 +01:00
|
|
|
std::cout << "-------------------- Unused parameters: --------------------\n";
|
|
|
|
|
param.displayUsage();
|
|
|
|
|
std::cout << "----------------------------------------------------------------" << std::endl;
|
2012-01-24 16:00:16 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Write parameters used for later reference.
|
2012-01-24 15:22:51 +01:00
|
|
|
if (output) {
|
2012-03-06 16:37:49 +01:00
|
|
|
param.writeParam(output_dir + "/spu_2p.param");
|
2012-01-24 15:22:51 +01:00
|
|
|
}
|
|
|
|
|
|
2012-01-23 13:06:25 +01:00
|
|
|
// Main simulation loop.
|
2012-02-09 23:15:14 +01:00
|
|
|
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();
|
2012-01-24 15:22:51 +01:00
|
|
|
std::cout << "\n\n================ Starting main simulation loop ===============" << std::endl;
|
2012-01-19 13:58:25 +01:00
|
|
|
for (int pstep = 0; pstep < num_psteps; ++pstep) {
|
2012-01-27 09:19:59 +01:00
|
|
|
std::cout << "\n\n--------------- Simulation step number " << pstep
|
2012-03-06 16:37:49 +01:00
|
|
|
<< " ---------------"
|
|
|
|
|
<< "\n Current time (days) " << Opm::unit::convert::to(current_time, Opm::unit::day)
|
|
|
|
|
<< "\n Current stepsize (days) " << Opm::unit::convert::to(stepsize, Opm::unit::day)
|
|
|
|
|
<< "\n Total time (days) " << Opm::unit::convert::to(total_time, Opm::unit::day)
|
|
|
|
|
<< "\n" << std::endl;
|
|
|
|
|
|
|
|
|
|
if (output) {
|
|
|
|
|
outputState(grid->c_grid(), state, pstep, output_dir);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Solve pressure.
|
|
|
|
|
if (use_gravity) {
|
|
|
|
|
computeTotalMobilityOmega(*props, allcells, state.saturation(), totmob, omega);
|
|
|
|
|
} else {
|
|
|
|
|
computeTotalMobility(*props, allcells, state.saturation(), totmob);
|
|
|
|
|
}
|
|
|
|
|
pressure_timer.start();
|
|
|
|
|
psolver.solve(totmob, omega, src, state.pressure(), state.faceflux());
|
|
|
|
|
pressure_timer.stop();
|
|
|
|
|
double pt = pressure_timer.secsSinceStart();
|
|
|
|
|
std::cout << "Pressure solver took: " << pt << " seconds." << std::endl;
|
|
|
|
|
ptime += pt;
|
|
|
|
|
|
|
|
|
|
// Solve transport
|
|
|
|
|
transport_timer.start();
|
|
|
|
|
if (use_reorder) {
|
|
|
|
|
// We must treat reorder_src here,
|
|
|
|
|
// if we are to handle anything but simple water
|
|
|
|
|
// injection, since it is expected to be
|
|
|
|
|
// equal to total outflow (if negative)
|
|
|
|
|
// and water inflow (if positive).
|
|
|
|
|
// Also, for anything but noflow boundaries,
|
|
|
|
|
// boundary flows must be accumulated into
|
|
|
|
|
// source term following the same convention.
|
|
|
|
|
Opm::toWaterSat(state.saturation(), reorder_sat);
|
|
|
|
|
reorder_model.solve(&state.faceflux()[0], &reorder_src[0], stepsize, &reorder_sat[0]);
|
|
|
|
|
Opm::toBothSat(reorder_sat, state.saturation());
|
|
|
|
|
if (use_segregation_split) {
|
|
|
|
|
if (use_column_solver) {
|
|
|
|
|
colsolver.solve(columns, stepsize, state.saturation());
|
|
|
|
|
} else {
|
|
|
|
|
std::vector<double> fluxes = state.faceflux();
|
|
|
|
|
std::fill(state.faceflux().begin(), state.faceflux().end(), 0.0);
|
|
|
|
|
tsolver.solve(*grid->c_grid(), tsrc, stepsize, ctrl, state, linsolve, rpt);
|
|
|
|
|
std::cout << rpt;
|
|
|
|
|
state.faceflux() = fluxes;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
} else {
|
|
|
|
|
tsolver.solve(*grid->c_grid(), tsrc, stepsize, ctrl, state, linsolve, rpt);
|
|
|
|
|
std::cout << rpt;
|
|
|
|
|
}
|
|
|
|
|
transport_timer.stop();
|
|
|
|
|
double tt = transport_timer.secsSinceStart();
|
|
|
|
|
std::cout << "Transport solver took: " << tt << " seconds." << std::endl;
|
|
|
|
|
ttime += tt;
|
|
|
|
|
|
|
|
|
|
current_time += stepsize;
|
2012-01-19 13:58:25 +01:00
|
|
|
}
|
2012-02-09 23:15:14 +01:00
|
|
|
total_timer.stop();
|
|
|
|
|
|
|
|
|
|
std::cout << "\n\n================ End of simulation ===============\n"
|
2012-03-06 16:37:49 +01:00
|
|
|
<< "Total time taken: " << total_timer.secsSinceStart()
|
|
|
|
|
<< "\n Pressure time: " << ptime
|
|
|
|
|
<< "\n Transport time: " << ttime << std::endl;
|
2012-01-19 14:12:53 +01:00
|
|
|
|
2012-01-20 15:51:29 +01:00
|
|
|
if (output) {
|
2012-03-06 16:37:49 +01:00
|
|
|
outputState(grid->c_grid(), state, num_psteps, output_dir);
|
2012-01-20 15:51:29 +01:00
|
|
|
}
|
2011-10-06 17:25:13 +02:00
|
|
|
|
2011-10-05 20:06:03 +02:00
|
|
|
destroy_transport_source(tsrc);
|
|
|
|
|
}
|