/*
Copyright 2014, 2015 Dr. Markus Blatt - HPC-Simulation-Software & Services.
Copyright 2014 Statoil AS
Copyright 2015 NTNU
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 .
*/
#include "config.h"
#include
namespace Opm
{
namespace AutoDiffGrid
{
Eigen::Array
faceCellsToEigen(const UnstructuredGrid& grid)
{
typedef Eigen::Array TwoColInt;
return Eigen::Map(grid.face_cells, grid.number_of_faces, 2);
}
Eigen::Array
cellCentroidsZToEigen(const UnstructuredGrid& grid)
{
return Eigen::Map >
(grid.cell_centroids, grid.number_of_cells, grid.dimensions).rightCols<1>();
}
/*
SparseTableView cell2Faces(const UnstructuredGrid& grid)
{
return SparseTableView(grid.cell_faces, grid.cell_facepos, numCells(grid));
}
*/
void extractInternalFaces(const UnstructuredGrid& grid,
Eigen::Array& internal_faces,
Eigen::Array& nbi)
{
typedef Eigen::Array OneColBool;
typedef Eigen::Array TwoColInt;
typedef Eigen::Array TwoColBool;
TwoColInt nb = faceCellsToEigen(grid);
// std::cout << "nb = \n" << nb << std::endl;
// Extracts the internal faces of the grid.
// These are stored in internal_faces.
TwoColBool nbib = nb >= 0;
OneColBool ifaces = nbib.rowwise().all();
const int num_internal = ifaces.cast().sum();
// std::cout << num_internal << " internal faces." << std::endl;
nbi.resize(num_internal, 2);
internal_faces.resize(num_internal);
int fi = 0;
int nf = numFaces(grid);
for (int f = 0; f < nf; ++f) {
if (ifaces[f]) {
internal_faces[fi] = f;
nbi.row(fi) = nb.row(f);
++fi;
}
}
}
} // end namespace AutoDiffGrid
#ifdef HAVE_DUNE_CORNERPOINT
namespace AutoDiffGrid
{
ADFaceCellTraits::Type
faceCellsToEigen(const Dune::CpGrid& grid)
{
return Dune::cpgrid::FaceCellsContainerProxy(&grid);
}
Eigen::Array
cellCentroidsZToEigen(const Dune::CpGrid& grid)
{
// Create an Eigen array of appropriate size
int rows=numCells(grid);
Eigen::Array array(rows);
// Fill it with the z coordinate of the cell centroids.
for (int i=0; i& internal_faces,
Eigen::Array& nbi)
{
// Extracts the internal faces of the grid.
// These are stored in internal_faces.
int nf=numFaces(grid);
int num_internal=0;
for(int f=0; f=0 && grid.faceCell(f, 1)>=0) {
internal_faces[fi] = f;
nbi(fi,0) = grid.faceCell(f, 0);
nbi(fi,1) = grid.faceCell(f, 1);
++fi;
}
}
}
} // end namespace AutoDiffGrid
#endif // HAVE_DUNE_CORNERPOINT
} // end namespace Opm