Fixed namespacing issue

This commit is contained in:
Kjetil Olsen Lye
2012-04-10 14:47:29 +02:00
parent d02e8e8651
commit 84b5041487
3 changed files with 32 additions and 8 deletions

View File

@@ -3,8 +3,11 @@
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/core/WellsManager.hpp>
#include <opm/core/GridManager.hpp>
#include <opm/core/pressure/IncompTpfa.hpp>
#include <opm/core/fluid/IncompPropertiesFromDeck.hpp>
#include <opm/core/linalg/LinearSolverUmfpack.hpp>
#include "opm/core/newwells.h"
#include "opm/core/grid.h"
int main(int argc, char** argv) {
@@ -12,7 +15,7 @@ int main(int argc, char** argv) {
using namespace Opm;
ParameterGroup parameters( argc, argv, false );
std::string file_name = parameters.getDefault<std::string>("inputdeck", "data.data");
// Read input file
EclipseGridParser parser(file_name);
std::cout << "Done!" << std::endl;
@@ -22,6 +25,18 @@ int main(int argc, char** argv) {
// Finally handle the wells
WellsManager wells(parser, *grid.c_grid(), NULL);
std::vector<int> global_cells(grid.c_grid()->global_cell, grid.c_grid()->global_cell + grid.c_grid()->number_of_cells);
std::cout << "ahoi" << std::endl;
double gravity[3] = {0.0, 0.0, parameters.getDefault<double>("gravity", 0.0)};
IncompPropertiesFromDeck incomp_properties(parser, global_cells);
std::cout << "there" << std::endl;
LinearSolverUmfpack umfpack_solver;
std::cout << "here" << std::endl;
IncompTpfa pressure_solver(*grid.c_grid(), incomp_properties.permeability(),
gravity, umfpack_solver, wells.c_wells());
return 0;
}

View File

@@ -38,14 +38,18 @@ namespace Opm
/// and N == g.number_of_cells.
/// \param[in] gravity Gravity vector. If nonzero, the array should
/// have D elements.
/// \param[in] wells The wells argument. Will be used in solution,
/// is ignored if NULL
IncompTpfa::IncompTpfa(const UnstructuredGrid& g,
const double* permeability,
const double* gravity,
const LinearSolverInterface& linsolver)
const LinearSolverInterface& linsolver,
const struct Wells* wells)
: grid_(g),
linsolver_(linsolver),
htrans_(g.cell_facepos[ g.number_of_cells ]),
trans_ (g.number_of_faces)
trans_ (g.number_of_faces),
wells_(wells)
{
UnstructuredGrid* gg = const_cast<UnstructuredGrid*>(&grid_);
tpfa_htrans_compute(gg, permeability, &htrans_[0]);
@@ -110,7 +114,7 @@ namespace Opm
}
}
ifs_tpfa_forces F = { NULL, NULL, NULL, NULL, NULL };
ifs_tpfa_forces F = { NULL, NULL, wells_, NULL, NULL };
if (! src.empty()) { F.src = &src[0]; }
F.bc = bcs;
@@ -172,7 +176,7 @@ namespace Opm
}
}
ifs_tpfa_forces F = { NULL, NULL, NULL, NULL, NULL };
ifs_tpfa_forces F = { NULL, NULL, wells_, NULL, NULL };
if (! src.empty()) { F.src = &src[0]; }
F.bc = bcs;

View File

@@ -25,6 +25,7 @@
struct UnstructuredGrid;
struct ifs_tpfa_data;
struct Wells;
struct FlowBoundaryConditions;
namespace Opm
@@ -47,10 +48,13 @@ namespace Opm
/// and N == g.number_of_cells.
/// \param[in] gravity Gravity vector. If nonzero, the array should
/// have D elements.
/// \param[in] wells The wells argument. Will be used in solution,
/// is ignored if NULL
IncompTpfa(const UnstructuredGrid& g,
const double* permeability,
const double* gravity,
const LinearSolverInterface& linsolver);
const LinearSolverInterface& linsolver,
const struct Wells* wells = 0);
/// Destructor.
~IncompTpfa();
@@ -113,7 +117,8 @@ namespace Opm
::std::vector<double> trans_ ;
::std::vector<double> gpress_;
::std::vector<double> gpress_omegaweighted_;
const struct Wells* wells_;
struct ifs_tpfa_data* h_;
};