mirror of
https://github.com/OPM/ResInsight.git
synced 2025-02-25 18:55:39 -06:00
#3927 Measurements : Add unit test for area computations
This commit is contained in:
@@ -20,15 +20,15 @@
|
||||
#include "gtest/gtest.h"
|
||||
|
||||
#include "cvfLibCore.h"
|
||||
#include "cvfLibViewing.h"
|
||||
#include "cvfLibRender.h"
|
||||
#include "cvfLibGeometry.h"
|
||||
#include "cvfLibRender.h"
|
||||
#include "cvfLibViewing.h"
|
||||
|
||||
#include "cvfArrayWrapperToEdit.h"
|
||||
#include "cvfArrayWrapperConst.h"
|
||||
#include "cvfArrayWrapperToEdit.h"
|
||||
|
||||
#include "cvfGeometryTools.h"
|
||||
#include "cvfBoundingBoxTree.h"
|
||||
#include "cvfGeometryTools.h"
|
||||
|
||||
#include <array>
|
||||
|
||||
@@ -86,10 +86,8 @@ void ControlVolume::calculateCubeFaceStatus(const cvf::Vec3dArray& nodeCoords, d
|
||||
}
|
||||
#endif
|
||||
|
||||
|
||||
template<typename NodeArrayType, typename NodeType, typename IndexType>
|
||||
NodeType quadNormal (ArrayWrapperConst<NodeArrayType, NodeType> nodeCoords,
|
||||
const IndexType cubeFaceIndices[4] )
|
||||
NodeType quadNormal(ArrayWrapperConst<NodeArrayType, NodeType> nodeCoords, const IndexType cubeFaceIndices[4])
|
||||
{
|
||||
return (nodeCoords[cubeFaceIndices[2]] - nodeCoords[cubeFaceIndices[0]]) ^
|
||||
(nodeCoords[cubeFaceIndices[3]] - nodeCoords[cubeFaceIndices[1]]);
|
||||
@@ -103,6 +101,7 @@ std::vector<cvf::Vec3d> createVertices()
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.resize(14, cvf::Vec3d::ZERO);
|
||||
|
||||
// clang-format off
|
||||
vxs[ 0]= cvf::Vec3d( 0 , 0 , 0 );
|
||||
vxs[ 1]= cvf::Vec3d( 1 , 0 , 0 );
|
||||
vxs[ 2]= cvf::Vec3d( 1 , 1 , 0 );
|
||||
@@ -117,6 +116,7 @@ std::vector<cvf::Vec3d> createVertices()
|
||||
vxs[11]= cvf::Vec3d( 1.2 , 0.6 , 0.0 );
|
||||
vxs[12]= cvf::Vec3d( 1.6 , 0.2 , 0.0 );
|
||||
vxs[13]= cvf::Vec3d( 0.8 ,-0.6 , 0.0 );
|
||||
// clang-format on
|
||||
|
||||
return vxs;
|
||||
}
|
||||
@@ -161,8 +161,7 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
edgeIntersectionStorage.setVertexCount(nodes.size());
|
||||
{
|
||||
std::vector<cvf::uint> polygon;
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(
|
||||
&polygon,
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
@@ -175,13 +174,11 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
EXPECT_TRUE(isOk);
|
||||
overlapPolygons.push_back(polygon);
|
||||
std::cout << polygon << std::endl;
|
||||
|
||||
}
|
||||
|
||||
{
|
||||
std::vector<cvf::uint> polygon;
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(
|
||||
&polygon,
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
@@ -194,13 +191,11 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
EXPECT_TRUE(isOk);
|
||||
overlapPolygons.push_back(polygon);
|
||||
std::cout << polygon << std::endl;
|
||||
|
||||
}
|
||||
|
||||
{
|
||||
std::vector<cvf::uint> polygon;
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(
|
||||
&polygon,
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
@@ -221,12 +216,7 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
|
||||
for (cvf::uint vxIdx = 0; vxIdx < nodes.size(); ++vxIdx)
|
||||
{
|
||||
GeometryTools::insertVertexInPolygon(
|
||||
&basePolygon,
|
||||
wrapArrayConst(&nodes),
|
||||
vxIdx,
|
||||
1e-6
|
||||
);
|
||||
GeometryTools::insertVertexInPolygon(&basePolygon, wrapArrayConst(&nodes), vxIdx, 1e-6);
|
||||
}
|
||||
|
||||
EXPECT_EQ((size_t)8, basePolygon.size());
|
||||
@@ -236,12 +226,7 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
{
|
||||
for (cvf::uint vxIdx = 0; vxIdx < nodes.size(); ++vxIdx)
|
||||
{
|
||||
GeometryTools::insertVertexInPolygon(
|
||||
&overlapPolygons[pIdx],
|
||||
wrapArrayConst(&nodes),
|
||||
vxIdx,
|
||||
1e-6
|
||||
);
|
||||
GeometryTools::insertVertexInPolygon(&overlapPolygons[pIdx], wrapArrayConst(&nodes), vxIdx, 1e-6);
|
||||
}
|
||||
|
||||
if (pIdx == 0)
|
||||
@@ -260,7 +245,6 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
std::cout << "Op" << pIdx << ":" << overlapPolygons[pIdx] << std::endl;
|
||||
}
|
||||
|
||||
|
||||
Vec3d normal = quadNormal(wrapArrayConst(&nodes), faces[0].data());
|
||||
std::vector<bool> faceOverlapPolygonWindingSameAsCubeFaceFlags;
|
||||
faceOverlapPolygonWindingSameAsCubeFaceFlags.resize(overlapPolygons.size(), true);
|
||||
@@ -275,15 +259,13 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
overlapPolygonPtrs.push_back(&(overlapPolygons[pIdx]));
|
||||
}
|
||||
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(
|
||||
wrapArrayConst(&nodes),
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(wrapArrayConst(&nodes),
|
||||
wrapArrayConst(&basePolygon),
|
||||
normal,
|
||||
overlapPolygonPtrs,
|
||||
faceOverlapPolygonWindingSameAsCubeFaceFlags,
|
||||
&freeFacePolygon,
|
||||
&hasHoles
|
||||
);
|
||||
&hasHoles);
|
||||
|
||||
EXPECT_EQ((size_t)4, freeFacePolygon.size());
|
||||
EXPECT_FALSE(hasHoles);
|
||||
@@ -300,24 +282,19 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
overlapPolygonPtrs.push_back(&(overlapPolygons[pIdx]));
|
||||
}
|
||||
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(
|
||||
wrapArrayConst(&nodes),
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(wrapArrayConst(&nodes),
|
||||
wrapArrayConst(&basePolygon),
|
||||
normal,
|
||||
overlapPolygonPtrs,
|
||||
faceOverlapPolygonWindingSameAsCubeFaceFlags,
|
||||
&freeFacePolygon,
|
||||
&hasHoles
|
||||
);
|
||||
&hasHoles);
|
||||
|
||||
EXPECT_EQ((size_t)9, freeFacePolygon.size());
|
||||
EXPECT_FALSE(hasHoles);
|
||||
|
||||
std::cout << "FF2: " << freeFacePolygon << std::endl;
|
||||
|
||||
}
|
||||
|
||||
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
@@ -326,12 +303,10 @@ TEST(CellFaceIntersectionTst, Intersection1)
|
||||
|
||||
TEST(CellFaceIntersectionTst, Intersection)
|
||||
{
|
||||
|
||||
std::vector<cvf::Vec3d> additionalVertices;
|
||||
cvf::Vec3dArray nodes;
|
||||
std::vector<size_t> polygon;
|
||||
|
||||
|
||||
cvf::Array<size_t> ids;
|
||||
size_t cv1CubeFaceIndices[4] = {0, 1, 2, 3};
|
||||
size_t cv2CubeFaceIndices[4] = {4, 5, 6, 7};
|
||||
@@ -352,9 +327,13 @@ TEST(CellFaceIntersectionTst, Intersection)
|
||||
nodes[6] = cvf::Vec3d(1, 1, 0);
|
||||
nodes[7] = cvf::Vec3d(0, 1, 0);
|
||||
|
||||
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon, &additionalVertices, &edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes), cv1CubeFaceIndices, cv2CubeFaceIndices, 1e-6);
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
cv1CubeFaceIndices,
|
||||
cv2CubeFaceIndices,
|
||||
1e-6);
|
||||
EXPECT_EQ((size_t)4, polygon.size());
|
||||
EXPECT_EQ((size_t)0, additionalVertices.size());
|
||||
EXPECT_TRUE(isOk);
|
||||
@@ -371,28 +350,27 @@ TEST(CellFaceIntersectionTst, Intersection)
|
||||
nodes[7] = cvf::Vec3d(-0.25, 0.5, 0);
|
||||
polygon.clear();
|
||||
|
||||
isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon, &additionalVertices, &edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes), cv1CubeFaceIndices, cv2CubeFaceIndices, 1e-6);
|
||||
isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
cv1CubeFaceIndices,
|
||||
cv2CubeFaceIndices,
|
||||
1e-6);
|
||||
EXPECT_EQ((size_t)8, polygon.size());
|
||||
EXPECT_EQ((size_t)8, additionalVertices.size());
|
||||
EXPECT_TRUE(isOk);
|
||||
|
||||
|
||||
|
||||
}
|
||||
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
TEST(CellFaceIntersectionTst, FreeFacePolygon)
|
||||
{
|
||||
|
||||
std::vector<cvf::Vec3d> additionalVertices;
|
||||
cvf::Vec3dArray nodes;
|
||||
std::vector<size_t> polygon;
|
||||
|
||||
|
||||
cvf::Array<size_t> ids;
|
||||
size_t cv1CubeFaceIndices[4] = {0, 1, 2, 3};
|
||||
size_t cv2CubeFaceIndices[4] = {4, 5, 6, 7};
|
||||
@@ -413,9 +391,13 @@ TEST(CellFaceIntersectionTst, FreeFacePolygon)
|
||||
nodes[6] = cvf::Vec3d(1, 1, 0);
|
||||
nodes[7] = cvf::Vec3d(0, 1, 0);
|
||||
|
||||
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon, &additionalVertices, &edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes), cv1CubeFaceIndices, cv2CubeFaceIndices, 1e-6);
|
||||
bool isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
cv1CubeFaceIndices,
|
||||
cv2CubeFaceIndices,
|
||||
1e-6);
|
||||
EXPECT_EQ((size_t)4, polygon.size());
|
||||
EXPECT_EQ((size_t)0, additionalVertices.size());
|
||||
EXPECT_TRUE(isOk);
|
||||
@@ -427,8 +409,7 @@ TEST(CellFaceIntersectionTst, FreeFacePolygon)
|
||||
|
||||
std::vector<size_t> partialFacePolygon;
|
||||
bool hasHoles = false;
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(
|
||||
wrapArrayConst(&nodes),
|
||||
GeometryTools::calculatePartiallyFreeCubeFacePolygon(wrapArrayConst(&nodes),
|
||||
wrapArrayConst(cv1CubeFaceIndices, 4),
|
||||
Vec3d(0, 0, 1),
|
||||
faceOverlapPolygons,
|
||||
@@ -448,12 +429,105 @@ TEST(CellFaceIntersectionTst, FreeFacePolygon)
|
||||
nodes[7] = cvf::Vec3d(-0.25, 0.5, 0);
|
||||
polygon.clear();
|
||||
|
||||
isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon, &additionalVertices, &edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes), cv1CubeFaceIndices, cv2CubeFaceIndices, 1e-6);
|
||||
isOk = GeometryTools::calculateOverlapPolygonOfTwoQuads(&polygon,
|
||||
&additionalVertices,
|
||||
&edgeIntersectionStorage,
|
||||
wrapArrayConst(&nodes),
|
||||
cv1CubeFaceIndices,
|
||||
cv2CubeFaceIndices,
|
||||
1e-6);
|
||||
EXPECT_EQ((size_t)8, polygon.size());
|
||||
EXPECT_EQ((size_t)8, additionalVertices.size());
|
||||
EXPECT_TRUE(isOk);
|
||||
|
||||
|
||||
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
///
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
TEST(CellFaceIntersectionTst, PolygonAreaNormal3D)
|
||||
{
|
||||
// Test special cases with zero area
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_TRUE(area == cvf::Vec3d::ZERO);
|
||||
}
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_TRUE(area == cvf::Vec3d::ZERO);
|
||||
}
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 1});
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_TRUE(area == cvf::Vec3d::ZERO);
|
||||
}
|
||||
|
||||
// Three points
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 1});
|
||||
vxs.push_back({0, 1, 1});
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_DOUBLE_EQ(-0.5, area.x());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.y());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.z());
|
||||
}
|
||||
|
||||
// four identical points
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 0});
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_TRUE(area == cvf::Vec3d::ZERO);
|
||||
}
|
||||
|
||||
// Square of four points
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 1});
|
||||
vxs.push_back({0, 1, 1});
|
||||
vxs.push_back({0, 1, 0});
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_DOUBLE_EQ(-1.0, area.x());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.y());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.z());
|
||||
}
|
||||
|
||||
// Square of four points + one point in center of square
|
||||
|
||||
{
|
||||
std::vector<cvf::Vec3d> vxs;
|
||||
vxs.push_back({0, 0, 0});
|
||||
vxs.push_back({0, 0, 1});
|
||||
vxs.push_back({0, 1, 1});
|
||||
vxs.push_back({0, 1, 0});
|
||||
|
||||
vxs.push_back({0, 0.5, 0.5}); // center of square
|
||||
|
||||
cvf::Vec3d area = GeometryTools::polygonAreaNormal3D(vxs);
|
||||
EXPECT_DOUBLE_EQ(-0.75, area.x());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.y());
|
||||
EXPECT_DOUBLE_EQ(0.0, area.z());
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user