[bugfix] fix build of ebos for Grid != CpGrid.

This commit is contained in:
Robert Kloefkorn 2017-06-22 15:55:01 +02:00
parent 18dd23ff44
commit 2aada96d8c
5 changed files with 79 additions and 19 deletions

View File

@ -184,17 +184,21 @@ namespace Ewoms
// no need to collect anything.
return;
}
typedef typename GridManager::GridView GridView;
typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView, Dune::MCMGElementLayout>
ElementMapper;
const CollectiveCommunication& comm = gridManager.grid().comm();
ElementMapper elemMapper(gridManager.gridView());
{
std::set< int > send, recv;
// the I/O rank receives from all other ranks
if( isIORank() )
{
typedef typename GridManager::EquilGrid::LeafGridView EquilGridView;
const EquilGridView equilGridView = gridManager.equilGrid().leafGridView() ;
typedef Dune::MultipleCodimMultipleGeomTypeMapper< EquilGridView, Dune::MCMGElementLayout>
EquilElementMapper;
EquilElementMapper equilElemMapper( equilGridView );
// the I/O rank needs a picture of the global grid, here we
// use equilGrid which represents a view on the global grid
const size_t globalSize = gridManager.equilGrid().leafGridView().size( 0 );
@ -206,9 +210,9 @@ namespace Ewoms
const auto& elemEndIt = gridManager.equilGrid().leafGridView().template end<0>();
for (; elemIt != elemEndIt; ++elemIt) {
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2,4)
int elemIdx = elemMapper.index(*elemIt );
int elemIdx = equilElemMapper.index(*elemIt );
#else
int elemIdx = elemMapper.map(*elemIt );
int elemIdx = equilElemMapper.map(*elemIt );
#endif
int cartElemIdx = gridManager.equilCartesianIndexMapper().cartesianIndex(elemIdx);
globalCartesianIndex_[elemIdx] = cartElemIdx;
@ -235,9 +239,16 @@ namespace Ewoms
IndexMapType distributedCartesianIndex;
distributedCartesianIndex.resize(gridSize, -1);
auto localView = gridManager.grid().leafGridView();
for( auto it = localView.template begin< 0 >(),
end = localView.template end< 0 >(); it != end; ++it )
typedef typename GridManager::GridView LocalGridView;
const LocalGridView localGridView = gridManager.gridView() ;
typedef Dune::MultipleCodimMultipleGeomTypeMapper< LocalGridView, Dune::MCMGElementLayout>
ElementMapper;
ElementMapper elemMapper( localGridView );
for( auto it = localGridView.template begin< 0, Dune::Interior_Partition >(),
end = localGridView.template end< 0, Dune::Interior_Partition >(); it != end; ++it )
{
const auto element = *it ;
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2,4)

View File

@ -161,7 +161,7 @@ public:
protected:
void createGrids_()
{
const auto& gridProps = this->eclState()->get3DProperties();
const auto& gridProps = this->eclState().get3DProperties();
const std::vector<double>& porv = gridProps.getDoubleGridProperty("PORV").getData();
// we use separate grid objects: one for the calculation of the initial condition
@ -173,7 +173,7 @@ protected:
// create the EQUIL grid
/////
equilGrid_ = new EquilGrid();
equilGrid_->processEclipseFormat(*this->eclState()->getInputGrid(),
equilGrid_->processEclipseFormat(this->eclState().getInputGrid(),
/*isPeriodic=*/false,
/*flipNormals=*/false,
/*clipZ=*/false,

View File

@ -144,10 +144,10 @@ public:
protected:
void createGrids_()
{
const auto& gridProps = this->eclState()->get3DProperties();
const auto& gridProps = this->eclState().get3DProperties();
const std::vector<double>& porv = gridProps.getDoubleGridProperty("PORV").getData();
grid_ = new Grid(*(this->deck()), porv);
grid_ = new Grid(this->deck(), porv);
cartesianIndexMapper_ = new CartesianIndexMapper(*grid_);
}

View File

@ -28,11 +28,14 @@
#ifndef EWOMS_ECL_PROBLEM_HH
#define EWOMS_ECL_PROBLEM_HH
//#define DISABLE_ALUGRID_SFC_ORDERING 1
//#define EBOS_USE_ALUGRID 1
// make sure that the EBOS_USE_ALUGRID macro. using the preprocessor for this is slightly
// hacky...
#if EBOS_USE_ALUGRID
//#define DISABLE_ALUGRID_SFC_ORDERING 1
#if !HAVE_DUNE_ALUGRID || !DUNE_VERSION_NEWER(DUNE_ALUGRID, 2,4)
#if !HAVE_DUNE_ALUGRID
#warning "ALUGrid was indicated to be used for the ECL black oil simulator, but this "
#warning "requires the presence of dune-alugrid >= 2.4. Falling back to Dune::CpGrid"
#undef EBOS_USE_ALUGRID

View File

@ -28,6 +28,7 @@
#ifndef EWOMS_ECL_TRANSMISSIBILITY_HH
#define EWOMS_ECL_TRANSMISSIBILITY_HH
#include <ewoms/common/propertysystem.hh>
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
@ -41,6 +42,8 @@
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/grid/CpGrid.hpp>
#include <array>
#include <vector>
#include <unordered_map>
@ -62,6 +65,7 @@ NEW_PROP_TAG(ElementMapper);
template <class TypeTag>
class EclTransmissibility
{
typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, GridManager) GridManager;
@ -96,6 +100,7 @@ public:
void finishInit()
{ update(); }
void update()
{
const auto& gridView = gridManager_.gridView();
@ -185,10 +190,14 @@ public:
unsigned insideFaceIdx = intersection.indexInInside();
unsigned outsideFaceIdx = intersection.indexInOutside();
int faceIdx = intersection.id();
DimVector faceCenterInside = gridManager_.grid().faceCenterEcl(insideElemIdx,insideFaceIdx);
DimVector faceCenterOutside = gridManager_.grid().faceCenterEcl(outsideElemIdx,outsideFaceIdx);
DimVector faceAreaNormal = gridManager_.grid().faceAreaNormalEcl(faceIdx);
DimVector faceCenterInside;
DimVector faceCenterOutside;
DimVector faceAreaNormal;
typename std::is_same< Grid, Dune::CpGrid> :: type isCpGrid;
computeFaceProperties( intersection, insideElemIdx, insideFaceIdx, outsideElemIdx, outsideFaceIdx,
faceCenterInside, faceCenterOutside, faceAreaNormal,
isCpGrid );
Scalar halfTrans1;
Scalar halfTrans2;
@ -266,6 +275,43 @@ public:
{ return trans_.at(isId_(elemIdx1, elemIdx2)); }
private:
template <class Intersection>
void computeFaceProperties( const Intersection& intersection,
const int insideElemIdx,
const int insideFaceIdx,
const int outsideElemIdx,
const int outsideFaceIdx,
DimVector& faceCenterInside,
DimVector& faceCenterOutside,
DimVector& faceAreaNormal,
std::false_type ) const
{
// default implementation for DUNE grids
const auto& geometry = intersection.geometry();
faceCenterInside = geometry.center();
faceCenterOutside = faceCenterInside;
faceAreaNormal = intersection.centerUnitOuterNormal();
faceAreaNormal *= geometry.volume();
}
template <class Intersection>
void computeFaceProperties( const Intersection& intersection,
const int insideElemIdx,
const int insideFaceIdx,
const int outsideElemIdx,
const int outsideFaceIdx,
DimVector& faceCenterInside,
DimVector& faceCenterOutside,
DimVector& faceAreaNormal,
std::true_type ) const
{
int faceIdx = intersection.id();
faceCenterInside = gridManager_.grid().faceCenterEcl(insideElemIdx,insideFaceIdx);
faceCenterOutside = gridManager_.grid().faceCenterEcl(outsideElemIdx,outsideFaceIdx);
faceAreaNormal = gridManager_.grid().faceAreaNormalEcl(faceIdx);
}
void extractPermeability_()
{
const auto& props = gridManager_.eclState().get3DProperties();