mirror of
https://github.com/OPM/ResInsight.git
synced 2025-02-25 18:55:39 -06:00
#5273 Allen Diagrams: Compute complete set of NNCs
Add class RigNncConnection Implement algorithm to compute the complete set of Nncs
This commit is contained in:
parent
360893817e
commit
fa791d0568
@ -177,8 +177,17 @@ void RivNNCGeometryGenerator::textureCoordinates( cvf::Vec2fArray*
|
||||
#pragma omp parallel for
|
||||
for ( int tIdx = 0; tIdx < static_cast<int>( m_triangleIndexToNNCIndex->size() ); tIdx++ )
|
||||
{
|
||||
double cellScalarValue = ( *nncResultVals )[( *m_triangleIndexToNNCIndex )[tIdx]];
|
||||
cvf::Vec2f texCoord = mapper->mapToTextureCoord( cellScalarValue );
|
||||
double cellScalarValue = HUGE_VAL;
|
||||
size_t resultIndex = ( *m_triangleIndexToNNCIndex )[tIdx];
|
||||
|
||||
// The nnc connections can have more connections than reported from Eclipse, clamp the result index to Eclipse Results
|
||||
|
||||
if ( resultIndex < nncResultVals->size() )
|
||||
{
|
||||
cellScalarValue = ( *nncResultVals )[resultIndex];
|
||||
}
|
||||
|
||||
cvf::Vec2f texCoord = mapper->mapToTextureCoord( cellScalarValue );
|
||||
if ( cellScalarValue == HUGE_VAL || cellScalarValue != cellScalarValue ) // a != a is true for NAN's
|
||||
{
|
||||
texCoord[1] = 1.0f;
|
||||
|
@ -76,6 +76,7 @@ ${CMAKE_CURRENT_LIST_DIR}/RigEquil.h
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigWbsParameter.h
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigEclipseAllenFaultsStatCalc.h
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigCellFaceGeometryTools.h
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigNncConnection.h
|
||||
)
|
||||
|
||||
|
||||
@ -149,6 +150,7 @@ ${CMAKE_CURRENT_LIST_DIR}/RigEquil.cpp
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigWbsParameter.cpp
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigEclipseAllenFaultsStatCalc.cpp
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigCellFaceGeometryTools.cpp
|
||||
${CMAKE_CURRENT_LIST_DIR}/RigNncConnection.cpp
|
||||
)
|
||||
|
||||
list(APPEND CODE_HEADER_FILES
|
||||
|
@ -20,9 +20,14 @@
|
||||
|
||||
#include "RigCell.h"
|
||||
#include "RigMainGrid.h"
|
||||
#include "RigNncConnection.h"
|
||||
|
||||
#include "cvfGeometryTools.h"
|
||||
|
||||
#include "cafAssert.h"
|
||||
|
||||
#include <QDebug>
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
@ -110,3 +115,227 @@ cvf::StructGridInterface::FaceType
|
||||
|
||||
return cvf::StructGridInterface::NO_FACE;
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
std::vector<RigConnection> RigCellFaceGeometryTools::computeOtherNncs( const RigMainGrid* mainGrid,
|
||||
const std::vector<RigConnection>& nativeConnections )
|
||||
{
|
||||
// Compute Non-Neighbor Connections (NNC) not reported by Eclipse. NNCs with zero transmissibility are not reported
|
||||
// by Eclipse. Use faults as basis for subset of cells to find NNC connection for. The imported connections from
|
||||
// Eclipse are located at the beginning of the connections vector.
|
||||
|
||||
std::vector<RigConnection> otherConnections;
|
||||
|
||||
class CellPair
|
||||
{
|
||||
public:
|
||||
CellPair( size_t globalIdx1, size_t globalIdx2 )
|
||||
{
|
||||
if ( globalIdx1 < globalIdx2 )
|
||||
{
|
||||
m_globalCellIdx1 = globalIdx1;
|
||||
m_globalCellIdx2 = globalIdx2;
|
||||
}
|
||||
else
|
||||
{
|
||||
m_globalCellIdx1 = globalIdx2;
|
||||
m_globalCellIdx2 = globalIdx1;
|
||||
}
|
||||
}
|
||||
|
||||
bool operator<( const CellPair& other ) const
|
||||
{
|
||||
if ( m_globalCellIdx1 != other.m_globalCellIdx1 )
|
||||
{
|
||||
return m_globalCellIdx1 < other.m_globalCellIdx1;
|
||||
}
|
||||
|
||||
return ( m_globalCellIdx2 < other.m_globalCellIdx2 );
|
||||
}
|
||||
|
||||
private:
|
||||
size_t m_globalCellIdx1;
|
||||
size_t m_globalCellIdx2;
|
||||
};
|
||||
|
||||
std::set<CellPair> nativeCellPairs;
|
||||
|
||||
for ( const auto& c : nativeConnections )
|
||||
{
|
||||
nativeCellPairs.emplace( CellPair( c.m_c1GlobIdx, c.m_c2GlobIdx ) );
|
||||
}
|
||||
|
||||
if ( nativeConnections.size() != nativeCellPairs.size() )
|
||||
{
|
||||
QString message = QString( "Nnc connection imported from Eclipse are not unique\nNNC count : %1\nUnique : %2" )
|
||||
.arg( nativeConnections.size() )
|
||||
.arg( nativeCellPairs.size() );
|
||||
qDebug() << message;
|
||||
}
|
||||
|
||||
std::set<CellPair> otherCellPairs;
|
||||
|
||||
const cvf::Collection<RigFault>& faults = mainGrid->faults();
|
||||
for ( size_t faultIdx = 0; faultIdx < faults.size(); faultIdx++ )
|
||||
{
|
||||
const RigFault* fault = faults.at( faultIdx );
|
||||
|
||||
const std::vector<RigFault::FaultFace>& faultFaces = fault->faultFaces();
|
||||
for ( const auto& f : faultFaces )
|
||||
{
|
||||
size_t sourceReservoirCellIndex = f.m_nativeReservoirCellIndex;
|
||||
cvf::StructGridInterface::FaceType sourceCellFace = f.m_nativeFace;
|
||||
|
||||
const std::vector<cvf::Vec3d>& mainGridNodes = mainGrid->nodes();
|
||||
|
||||
cvf::BoundingBox bb;
|
||||
std::array<size_t, 4> sourceFaceIndices;
|
||||
mainGrid->globalCellArray()[sourceReservoirCellIndex].faceIndices( sourceCellFace, &sourceFaceIndices );
|
||||
|
||||
bb.add( mainGridNodes[sourceFaceIndices[0]] );
|
||||
bb.add( mainGridNodes[sourceFaceIndices[1]] );
|
||||
bb.add( mainGridNodes[sourceFaceIndices[2]] );
|
||||
bb.add( mainGridNodes[sourceFaceIndices[3]] );
|
||||
|
||||
std::vector<size_t> closeCells;
|
||||
mainGrid->findIntersectingCells( bb, &closeCells );
|
||||
|
||||
cvf::StructGridInterface::FaceType candidateFace = cvf::StructGridInterface::oppositeFace( sourceCellFace );
|
||||
|
||||
size_t neighborCellIndex = std::numeric_limits<size_t>::max();
|
||||
size_t ni = std::numeric_limits<size_t>::max();
|
||||
size_t nj = std::numeric_limits<size_t>::max();
|
||||
size_t nk = std::numeric_limits<size_t>::max();
|
||||
|
||||
{
|
||||
size_t i;
|
||||
size_t j;
|
||||
size_t k;
|
||||
mainGrid->ijkFromCellIndex( sourceReservoirCellIndex, &i, &j, &k );
|
||||
|
||||
mainGrid->neighborIJKAtCellFace( i, j, k, sourceCellFace, &ni, &nj, &nk );
|
||||
|
||||
if ( mainGrid->isCellValid( ni, nj, nk ) )
|
||||
{
|
||||
neighborCellIndex = mainGrid->cellIndexFromIJK( ni, nj, nk );
|
||||
}
|
||||
}
|
||||
|
||||
for ( size_t candidateCellIndex : closeCells )
|
||||
{
|
||||
if ( candidateCellIndex == sourceReservoirCellIndex )
|
||||
{
|
||||
// Exclude cellIndex for source cell
|
||||
continue;
|
||||
}
|
||||
|
||||
if ( candidateCellIndex == neighborCellIndex )
|
||||
{
|
||||
// Exclude direct neighbor
|
||||
continue;
|
||||
}
|
||||
|
||||
if ( neighborCellIndex != std::numeric_limits<size_t>::max() )
|
||||
{
|
||||
// Find target IJK index based on source cell and cell face
|
||||
// Exclude cells not matching destination target index
|
||||
|
||||
size_t ci = std::numeric_limits<size_t>::max();
|
||||
size_t cj = std::numeric_limits<size_t>::max();
|
||||
size_t ck = std::numeric_limits<size_t>::max();
|
||||
mainGrid->ijkFromCellIndex( candidateCellIndex, &ci, &cj, &ck );
|
||||
|
||||
if ( sourceCellFace == cvf::StructGridInterface::POS_I ||
|
||||
sourceCellFace == cvf::StructGridInterface::NEG_I )
|
||||
{
|
||||
if ( ni != ci )
|
||||
{
|
||||
continue;
|
||||
}
|
||||
}
|
||||
else if ( sourceCellFace == cvf::StructGridInterface::POS_J ||
|
||||
sourceCellFace == cvf::StructGridInterface::NEG_J )
|
||||
{
|
||||
if ( nj != cj )
|
||||
{
|
||||
continue;
|
||||
}
|
||||
}
|
||||
else if ( sourceCellFace == cvf::StructGridInterface::POS_K ||
|
||||
sourceCellFace == cvf::StructGridInterface::NEG_K )
|
||||
{
|
||||
if ( nk != ck )
|
||||
{
|
||||
continue;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
CellPair candidate( sourceReservoirCellIndex, candidateCellIndex );
|
||||
|
||||
if ( nativeCellPairs.count( candidate ) > 0 )
|
||||
{
|
||||
continue;
|
||||
}
|
||||
|
||||
if ( otherCellPairs.count( candidate ) > 0 )
|
||||
{
|
||||
continue;
|
||||
}
|
||||
|
||||
std::vector<size_t> polygon;
|
||||
std::vector<cvf::Vec3d> intersections;
|
||||
|
||||
std::array<size_t, 4> candidateFaceIndices;
|
||||
mainGrid->globalCellArray()[candidateCellIndex].faceIndices( candidateFace, &candidateFaceIndices );
|
||||
|
||||
bool foundOverlap =
|
||||
cvf::GeometryTools::calculateOverlapPolygonOfTwoQuads( &polygon,
|
||||
&intersections,
|
||||
(cvf::EdgeIntersectStorage<size_t>*)nullptr,
|
||||
cvf::wrapArrayConst( &mainGridNodes ),
|
||||
sourceFaceIndices.data(),
|
||||
candidateFaceIndices.data(),
|
||||
1e-6 );
|
||||
|
||||
if ( foundOverlap )
|
||||
{
|
||||
otherCellPairs.emplace( candidate );
|
||||
|
||||
RigConnection conn;
|
||||
conn.m_c1GlobIdx = sourceReservoirCellIndex;
|
||||
conn.m_c1Face = sourceCellFace;
|
||||
conn.m_c2GlobIdx = candidateCellIndex;
|
||||
|
||||
conn.m_polygon = RigCellFaceGeometryTools::extractPolygon( mainGridNodes, polygon, intersections );
|
||||
|
||||
otherConnections.emplace_back( conn );
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
return otherConnections;
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
std::vector<cvf::Vec3d> RigCellFaceGeometryTools::extractPolygon( const std::vector<cvf::Vec3d>& nativeNodes,
|
||||
const std::vector<size_t>& connectionPolygon,
|
||||
const std::vector<cvf::Vec3d>& connectionIntersections )
|
||||
{
|
||||
std::vector<cvf::Vec3d> allPolygonNodes;
|
||||
|
||||
for ( size_t polygonIndex : connectionPolygon )
|
||||
{
|
||||
if ( polygonIndex < nativeNodes.size() )
|
||||
allPolygonNodes.push_back( nativeNodes[polygonIndex] );
|
||||
else
|
||||
allPolygonNodes.push_back( connectionIntersections[polygonIndex - nativeNodes.size()] );
|
||||
}
|
||||
|
||||
return allPolygonNodes;
|
||||
}
|
||||
|
@ -26,6 +26,7 @@
|
||||
|
||||
class RigCell;
|
||||
class RigMainGrid;
|
||||
class RigConnection;
|
||||
|
||||
//==================================================================================================
|
||||
///
|
||||
@ -38,4 +39,11 @@ public:
|
||||
const RigMainGrid& mainGrid,
|
||||
std::vector<size_t>* connectionPolygon,
|
||||
std::vector<cvf::Vec3d>* connectionIntersections );
|
||||
|
||||
static std::vector<RigConnection> computeOtherNncs( const RigMainGrid* mainGrid,
|
||||
const std::vector<RigConnection>& nativeConnections );
|
||||
|
||||
static std::vector<cvf::Vec3d> extractPolygon( const std::vector<cvf::Vec3d>& nativeNodes,
|
||||
const std::vector<size_t>& connectionPolygon,
|
||||
const std::vector<cvf::Vec3d>& connectionIntersections );
|
||||
};
|
||||
|
@ -556,6 +556,8 @@ void RigMainGrid::calculateFaults( const RigActiveCellInfo* activeCellInfo )
|
||||
}
|
||||
}
|
||||
|
||||
this->nncData()->computeNncsFromFaults( this );
|
||||
|
||||
distributeNNCsToFaults();
|
||||
}
|
||||
|
||||
|
@ -57,14 +57,10 @@ void RigNNCData::processConnections( const RigMainGrid& mainGrid )
|
||||
// Found an overlap polygon. Store data about connection
|
||||
|
||||
m_connections[cnIdx].m_c1Face = connectionFace;
|
||||
for ( size_t pIdx = 0; pIdx < connectionPolygon.size(); ++pIdx )
|
||||
{
|
||||
if ( connectionPolygon[pIdx] < mainGrid.nodes().size() )
|
||||
m_connections[cnIdx].m_polygon.push_back( mainGrid.nodes()[connectionPolygon[pIdx]] );
|
||||
else
|
||||
m_connections[cnIdx].m_polygon.push_back(
|
||||
connectionIntersections[connectionPolygon[pIdx] - mainGrid.nodes().size()] );
|
||||
}
|
||||
|
||||
m_connections[cnIdx].m_polygon = RigCellFaceGeometryTools::extractPolygon( mainGrid.nodes(),
|
||||
connectionPolygon,
|
||||
connectionIntersections );
|
||||
|
||||
// Add to search map, possibly not needed
|
||||
// m_cellIdxToFaceToConnectionIdxMap[m_connections[cnIdx].m_c1GlobIdx][connectionFace].push_back(cnIdx);
|
||||
@ -78,6 +74,18 @@ void RigNNCData::processConnections( const RigMainGrid& mainGrid )
|
||||
}
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
void RigNNCData::computeNncsFromFaults( const RigMainGrid* mainGrid )
|
||||
{
|
||||
m_nativeConnectionCount = m_connections.size();
|
||||
|
||||
std::vector<RigConnection> otherConnections = RigCellFaceGeometryTools::computeOtherNncs( mainGrid, m_connections );
|
||||
|
||||
m_connections.insert( m_connections.end(), otherConnections.begin(), otherConnections.end() );
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
|
@ -21,6 +21,8 @@
|
||||
|
||||
#include "RiaNncDefines.h"
|
||||
|
||||
#include "RigNncConnection.h"
|
||||
|
||||
#include "cvfObject.h"
|
||||
#include "cvfStructGrid.h"
|
||||
#include "cvfVector3.h"
|
||||
@ -33,28 +35,6 @@ class RigMainGrid;
|
||||
class RigCell;
|
||||
class RigEclipseResultAddress;
|
||||
|
||||
class RigConnection
|
||||
{
|
||||
public:
|
||||
RigConnection()
|
||||
: m_c1GlobIdx( cvf::UNDEFINED_SIZE_T )
|
||||
, m_c1Face( cvf::StructGridInterface::NO_FACE )
|
||||
, m_c2GlobIdx( cvf::UNDEFINED_SIZE_T )
|
||||
{
|
||||
}
|
||||
|
||||
bool hasCommonArea() const
|
||||
{
|
||||
return m_polygon.size() > 0;
|
||||
}
|
||||
|
||||
size_t m_c1GlobIdx;
|
||||
cvf::StructGridInterface::FaceType m_c1Face;
|
||||
size_t m_c2GlobIdx;
|
||||
|
||||
std::vector<cvf::Vec3d> m_polygon;
|
||||
};
|
||||
|
||||
class RigNNCData : public cvf::Object
|
||||
{
|
||||
public:
|
||||
@ -68,6 +48,7 @@ public:
|
||||
RigNNCData();
|
||||
|
||||
void processConnections( const RigMainGrid& mainGrid );
|
||||
void computeNncsFromFaults( const RigMainGrid* mainGrid );
|
||||
|
||||
void setConnections( std::vector<RigConnection>& connections );
|
||||
|
||||
@ -108,6 +89,7 @@ private:
|
||||
|
||||
private:
|
||||
std::vector<RigConnection> m_connections;
|
||||
size_t m_nativeConnectionCount;
|
||||
std::map<QString, std::vector<std::vector<double>>> m_connectionResults;
|
||||
std::map<RigEclipseResultAddress, QString> m_resultAddrToNNCDataType;
|
||||
};
|
||||
|
39
ApplicationCode/ReservoirDataModel/RigNncConnection.cpp
Normal file
39
ApplicationCode/ReservoirDataModel/RigNncConnection.cpp
Normal file
@ -0,0 +1,39 @@
|
||||
/////////////////////////////////////////////////////////////////////////////////
|
||||
//
|
||||
// Copyright (C) 2020 Equinor ASA
|
||||
//
|
||||
// ResInsight 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.
|
||||
//
|
||||
// ResInsight 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 at <http://www.gnu.org/licenses/gpl.html>
|
||||
// for more details.
|
||||
//
|
||||
/////////////////////////////////////////////////////////////////////////////////
|
||||
|
||||
#include "RigNncConnection.h"
|
||||
|
||||
#include "cvfMath.h"
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
RigConnection::RigConnection()
|
||||
: m_c1GlobIdx( cvf::UNDEFINED_SIZE_T )
|
||||
, m_c1Face( cvf::StructGridInterface::NO_FACE )
|
||||
, m_c2GlobIdx( cvf::UNDEFINED_SIZE_T )
|
||||
{
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
bool RigConnection::hasCommonArea() const
|
||||
{
|
||||
return m_polygon.size() > 0;
|
||||
}
|
41
ApplicationCode/ReservoirDataModel/RigNncConnection.h
Normal file
41
ApplicationCode/ReservoirDataModel/RigNncConnection.h
Normal file
@ -0,0 +1,41 @@
|
||||
/////////////////////////////////////////////////////////////////////////////////
|
||||
//
|
||||
// Copyright (C) 2020 Equinor ASA
|
||||
//
|
||||
// ResInsight 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.
|
||||
//
|
||||
// ResInsight 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 at <http://www.gnu.org/licenses/gpl.html>
|
||||
// for more details.
|
||||
//
|
||||
/////////////////////////////////////////////////////////////////////////////////
|
||||
|
||||
#pragma once
|
||||
|
||||
#include "cvfStructGrid.h"
|
||||
#include "cvfVector3.h"
|
||||
|
||||
#include <vector>
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
class RigConnection
|
||||
{
|
||||
public:
|
||||
RigConnection();
|
||||
|
||||
bool hasCommonArea() const;
|
||||
|
||||
size_t m_c1GlobIdx;
|
||||
cvf::StructGridInterface::FaceType m_c1Face;
|
||||
size_t m_c2GlobIdx;
|
||||
|
||||
std::vector<cvf::Vec3d> m_polygon;
|
||||
};
|
Loading…
Reference in New Issue
Block a user