3D UGrid Tutorial
Introduction
The purpose of this tutorial is to provide explanation on how to use xmsgrid to create three dimensional unstructured grids, or UGrids. A UGrid has points and cells defined using those points. There are many different kinds of cells available, but this tutorial will focus on the 3-dimensional cells, namely: tetrahedron, voxel, hexahedron, wedge, pyramid, and polyhedron.
Example - Defining Ugrid Cells
Supported 3D grid cells include: tetrahedron (10), voxel (11), hexahedron (12), wedge (13), pyramid (14), and polyhedron (42). A cell is defined with the number declaration of the shape (10, 11, 12, 13, 14, and 42, respectively, as defined by the enumeration xms::XmUGridCellType), then the number of points, followed by the point indices. The cell definitions mirror VTK cell definitions which are available on page 9 of VTK File Formats for VTK version 4.2 at https://www.vtk.org/wp-content/uploads/2015/04/file-formats.pdf.
A tetrahedron (10) has 4 points. See the illustration in the VTK file Format pdf referenced above on page 9 for a VTK_TETRA for point order. A cellstream example for a tetrahedron is: 10, 4, 0, 1, 2, 3.
A voxel (11) has 8 orthogonially defined points. See the illustration in the VTK file Format pdf referenced above on page 9 for a VTK_VOXEL for point order. A cellstream example for a voxel is: 11, 8, 0, 1, 2, 3, 4, 5, 6, 7.
A hexahedron (12) has 8 points. See the illustration in the VTK file Format pdf referenced above on page 9 for a VTK_HEXAHEDRON for point order. A cellstream example for a hexahedron is: 12, 8, 0, 1, 2, 3, 4, 5, 6, 7.
A wedge (13) has 6 points. See the illustration in the VTK file Format pdf referenced above on page 9 for a VTK_WEDGE for point order. A cellstream example for a wedge is: 13, 6, 0, 1, 2, 3, 4, 5.
A pyramid (14) has 5 points. See the illustration in the VTK file Format pdf referenced above on page 9 for a VTK_PYRAMID for point order. A cellstream example for a pyramid is: 14, 5, 0, 1, 2, 3, 4.
A polyhedron (42) does not a have a defined number of points. A polyhedron cellstream has the following format: The cell type, number of faces, and then repeated for each face: the number of points in the face, followed by the points in the face declared in a counter-clockwise (ccw) direction (as viewed from outside the solid). A cellstream example for a polyhedron representation of a cube defined with points (10, 11, 12, 13) on the top face and points (14, 15, 16, 17) on the bottom face as viewed from above with points ordered clockwise from the lower left as viewed from above) is: 42, 6, 4, 10, 13, 12, 11, 4, 10, 11, 15, 14, 4, 11, 12, 16, 15, 4, 12, 13, 17, 16, 4, 13, 10, 14, 17, 4, 14, 15, 16, 17. (42 is the type number, 6 is the number of faces, 4 is the number of points in the first face, those points are: 10, 13, 12, 11 in ccw order from the outside, 4 is the number of points on the next face, 10, 11, 15, 14 are the point indices of that face in ccw order, and so on.)
The testing code for this example is in TEST_XmUGrid3DLinear.
VecPt3d points = {{0, 0, 0}, {10, 0, 0}, {20, 0, 0}, {30, 0, 0}, {40, 0, 0},
{0, 10, 0}, {10, 10, 0}, {20, 10, 0}, {30, 10, 0}, {40, 10, 0},
{0, 20, 0}, {10, 20, 0}, {20, 20, 0}, {30, 20, 0}, {40, 20, 0},
{0, 0, 10}, {10, 0, 10}, {20, 0, 10}, {30, 0, 10}, {40, 0, 10},
{0, 10, 10}, {10, 10, 10}, {20, 10, 10}, {30, 10, 10}, {40, 10, 10},
{0, 20, 10}, {10, 20, 10}, {20, 20, 10}, {30, 20, 10}, {40, 20, 10}};
std::vector<int> cells = {(int)XMU_TETRA, 4, 0, 1, 5, 15,
(int)XMU_VOXEL, 8, 1, 2, 6, 7, 16, 17, 21, 22,
(int)XMU_HEXAHEDRON, 8, 2, 3, 8, 7, 17, 18, 23, 22,
(int)XMU_POLYHEDRON, 6,
4, 9, 8, 13, 14,
4, 8, 9, 24, 23,
4, 9, 14, 29, 24,
4, 14, 13, 28, 29,
4, 8, 13, 28, 23,
4, 23, 24, 29, 28,
(int)XMU_WEDGE, 6, 3, 4, 18, 8, 9, 23,
(int)XMU_PYRAMID, 5, 5, 6, 11, 10, 20};
Example - Creating a New 3D UGrid With Data
This example shows how to create a new 3D UGrid with the overloaded New function which directly sets the data as it constructs the grid, using the static function xms::XmUGrid::New(const VecPt3d& a_points, const VecInt& a_cellStream). Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::testGetSetPoint.
VecPt3d points = {{0, 10, 0}, {10, 10, 10}, {20, 10, 0}, {0, 0, 0}, {10, 0, 0},
{20, 0, 0}, {0, -10, 0}, {10, -10, 0}, {20, -10, 0}};
VecInt cellstream = {9, 4, 0, 3, 4, 1, 9, 4, 1, 4, 5, 2, 9, 4, 3, 6, 7, 4, 9, 4, 4, 7, 8, 5};
std::shared_ptr<XmUGrid> ugrid = XmUGrid::New(points, cellstream);
Example - Creating a New 3D UGrid
This example shows how to create a new 3D UGrid. XmUGrid is an abstract class that cannot be instantiated. The static function xms::XmUGrid::New is the only method to obtain a Boost Shared Pointer to an instance of XmUGrid. Points and Cellstream may be passed to New to initialize the UGrid with data but this is not required. Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::testUGridStreams.
Example - Setting the UGrid Points
This example shows how to set all of the UGrid points using xms::XmUGrid::SetLocations. The Ugrid points cannot be added or removed individually, though they can be edited individually. The SetPoints function takes one argument, a vector of 3D points. It is reccomended that each point is unique to avoid unexpected behavior. Functionality is shared between 2d and 3d UGrids. The testing code for this example is the same as used for creating a new UGrid, XmUGridUnitTests::testUGridStreams.
Example - Setting the UGrid Cell Stream
This example shows how to set the entire UGrid Cellstream using xms::XmUGrid::SetCellstream. Cellstreams are formatted as a stream of integers starting with the cell type as described by the enumeration XmUGridCellType, then the number of points in the cell, followed by a series of indices to points. The SetCellstream function takes one argument, a vector of integers formatted as previously described. Functionality is shared between 2d and 3d UGrids. The testing code for this example is the same as used for creating a new UGrid, XmUGridUnitTests::testUGridStreams.
The xms::XmUGrid::IsValidCellstream function checks whether all cell types and point counts within the cellstream match up, but does not check for orientation or valid point positions. This function should not be relied on to catch all errors, but can be a basic check before setting the cellstream. The function takes one arguement, a cellstream containing one or more cells. Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::SetCellstream.
{
std::shared_ptr<XmUGrid> emptyUGrid = XmUGrid::New();
VecPt3d points = emptyUGrid->GetLocations();
TS_ASSERT_EQUALS(0, points.size());
VecInt cellstream = emptyUGrid->GetCellstream();
TS_ASSERT(XmUGrid::IsValidCellstream(cellstream));
TS_ASSERT_EQUALS(0, cellstream.size());
std::shared_ptr<XmUGrid> ugrid = XmUGrid::New();
points = {{0, 10, 0}, {10, 10, 0}, {20, 10, 0}, {0, 0, 0}, {10, 0, 0},
{20, 0, 0}, {0, -10, 0}, {10, -10, 0}, {20, -10, 0}};
cellstream = {9, 4, 0, 3, 4, 1, 9, 4, 1, 4, 5, 2, 9, 4, 3, 6, 7, 4, 9, 4, 4, 7, 8, 5};
TS_ASSERT(!ugrid->GetModified());
ugrid->SetLocations(points);
VecPt3d pointsOut = ugrid->GetLocations();
TS_ASSERT_EQUALS(points, pointsOut);
TS_ASSERT(ugrid->GetModified());
ugrid->SetUnmodified();
TS_ASSERT(!ugrid->GetModified());
TS_ASSERT(ugrid->SetCellstream(cellstream));
VecInt cellstreamOut = ugrid->GetCellstream();
TS_ASSERT_EQUALS(cellstream, cellstreamOut);
TS_ASSERT(ugrid->GetModified());
cellstream = {-1};
TS_ASSERT(!XmUGrid::IsValidCellstream(cellstream));
cellstream = {9, 4, 0, 3, 4};
TS_ASSERT(!XmUGrid::IsValidCellstream(cellstream));
cellstream = {9, 3, 0, 3, 4, 1};
TS_ASSERT(!XmUGrid::IsValidCellstream(cellstream));
std::shared_ptr<XmUGrid> ugridBadOrder = XmUGrid::New();
TS_ASSERT(!ugridBadOrder->SetCellstream(cellstream));
}
Example - Get Number Of Points
This example shows how to return the number of points contained in a UGrid. The xms::XmUGrid::GetPointCount function returns the number of points in the UGrid. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get Point Locations
This example shows how to get all points contained within the XmUGrid. The xms::XmUGrid::GetLocations function returns a vector of Pt3d's. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get Location of a Point
This example shows how to get a specific point given a point index. The xms::XmUGrid::GetPointLocation function returns the Pt3d of the point at the specified index. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Set Point Location
This example shows how to set a specific point given a point index. The xms::XmUGrid::SetLocation function takes a point index and a Pt3d as arguements and returns whether the operation was succesful. The function may fail and return false if the change would cause any edges to cross. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get Locations of an Array of Points
This example shows how to convert a vector of point indices into a vector of Pt3d's. The xms::XmUGrid::GetPointsLocations function takes a vector of Point Indices as an arguement and returns a vector of corresponding Pt3d's as a result. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get GetExtents of UGrid
This example shows how to get the extents of a UGrid. The xms::XmUGrid::GetExtents function takes two 3D points as arguements describing the minimum and maximum extent that the points of the UGrid cover. These arguments will be set by the funtion. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get Cells Associated with a Point
This example shows how to the cells associated with a single point. The xms::XmUGrid::GetPointCells function takes one point index as an arguement and returns a vector of point indices. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
Example - Get the Cells that Share the Same Point or Points
This example shows how to the cells that share the same point or points. The xms::XmUGrid::GetPointsAdjacentCells function takes a vector of point indices as an arguement and returns a vector of cell indices. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testPointFunctions.
{
TS_ASSERT_EQUALS(ugrid->GetPointCount(), 9);
VecPt3d points = {{0, 10, 0}, {10, 10, 0}, {20, 10, 0}, {0, 0, 0}, {10, 0, 0},
{20, 0, 0}, {0, -10, 0}, {10, -10, 0}, {20, -10, 0}};
TS_ASSERT_EQUALS(ugrid->GetLocations(), points);
for (int i = 0; i < points.size(); ++i)
{
TS_ASSERT_EQUALS(points[i], ugrid->GetPointLocation(i));
}
Pt3d invalid = {100, 100, 100};
TS_ASSERT(!ugrid->IsValidPointChange(4, invalid));
TS_ASSERT(ugrid->SetPointLocation(4, valid));
TS_ASSERT_EQUALS(valid, ugrid->GetPointLocation(4));
VecInt pointIndices = {0, 3, 6};
VecPt3d expectedPoints = {{0, 10, 0}, {0, 0, 0}, {0, -10, 0}};
pointsPt3d = ugrid->GetPointsLocations(pointIndices);
TS_ASSERT_EQUALS(expectedPoints, pointsPt3d);
Pt3d min, max, expectedMin = {0, -10, 0}, expectedMax = {20, 10, 5};
ugrid->GetExtents(min, max);
TS_ASSERT_EQUALS(expectedMin, min);
TS_ASSERT_EQUALS(expectedMax, max);
VecInt expectedCells = {0, 2};
VecInt cellsAssociated = ugrid->GetPointAdjacentCells(3);
TS_ASSERT_EQUALS(expectedCells, cellsAssociated);
VecInt expectedCommonCells = {0, 1};
VecInt commonCells = ugrid->GetPointsAdjacentCells(pointIndex);
TS_ASSERT_EQUALS(expectedCells, cellsAssociated);
}
Example - Get the Points of a Cell
This example shows how to get the points of a cell. The xms::XmUGrid::GetCellPoints function takes one cell index and returns a vector of point indices. There is also a The xms::XmUGrid::GetCellLocations function that returns the locations (Pt3d) of those points. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellFunctions.
Example - Get the Type of a Cell
This example shows how to get the type of a cell. The xms::XmUGrid::GetCellType function takes one cell index and returns the cell type. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellFunctions.
Example - Get the Number of the Cells of each Dimension in a UGrid
This example shows how to get the count of the dimensions of cells used in a UGrid. The xms::XmUGrid::GetDimensionCounts function returns a vector with the number of cells with zero dimensions in the zero index, number of cells with one dimension in the first index, number of cells with two dimensions in the second index, and the number of cells with three dimensions in the third index. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellFunctions.
Example - Get the Dimension of a Cell
This example shows how to get the dimension of a cell. The xms::XmUGrid::GetCellDimension function takes one cell index and returns the cell dimension as an integer value. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellFunctions.
{
VecInt pointsOfCell = ugrid->GetCellPoints(0);
VecInt expectedPoints = {0, 3, 4, 1};
TS_ASSERT_EQUALS(expectedPoints, pointsOfCell);
ugrid->GetCellPoints(0, pointsOfCell);
TS_ASSERT_EQUALS(expectedPoints, pointsOfCell);
ugrid->GetCellLocations(0, locations);
VecPt3d expectedLocations = {{0, 10, 0}, {0, 0, 0}, {10, 0, 0}, {10, 10, 0}};
TS_ASSERT_EQUALS(expectedLocations, locations);
ugrid->GetCellExtents(0, min, max);
TS_ASSERT_EQUALS(
Pt3d(0, 0, 0), min);
TS_ASSERT_EQUALS(
Pt3d(10, 10, 0), max);
TS_ASSERT_EQUALS(XMU_QUAD, ugrid->GetCellType(0));
VecInt expectedDimensions = {0, 0, 4, 0};
TS_ASSERT_EQUALS(expectedDimensions, ugrid->GetDimensionCounts());
TS_ASSERT_EQUALS(2, ugrid->GetCellDimension(0));
TS_ASSERT(ugrid->GetCellCentroid(0, centroid));
TS_ASSERT_EQUALS(
Pt3d(5, 5, 0), centroid);
}
Example - Get the Cellstream of the UGrid
This example shows how to get the entire cellstream of the UGrid. The xms::XmUGrid::GetCellStream function returns a vector of integers that is the cellstream. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellstreamFunctions.
Example - Get a Single Cellstream for One Cell
This example shows how to get the cellstream for just one cell. The xms::XmUGrid::GetCellCellstream function takes a cell index and a vector of integers that is the cellstream for the specified cell passed in by reference. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellstreamFunctions.
{
VecInt cellstream = ugrid->GetCellstream();
VecInt expectedCellstream = {XMU_QUAD, 4, 0, 3, 4, 1, XMU_QUAD, 4, 1, 4, 5, 2,
XMU_QUAD, 4, 3, 6, 7, 4, XMU_QUAD, 4, 4, 7, 8, 5};
TS_ASSERT_EQUALS(expectedCellstream, cellstream);
ugrid->GetCellCellstream(0, cellstream);
expectedCellstream = {XMU_QUAD, 4, 0, 3, 4, 1};
TS_ASSERT_EQUALS(expectedCellstream, cellstream);
}
Example - Get the Cells Adjacent to a Given Cell
This example shows how to get the cells that also use any of the points of a given cell. Cells are adjacent if they share at least one point. The xms::XmUGrid::GetCellAdjacentCells function takes a cell index and returns a vector of cell indices. Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::testGetCellAdjacentCells.
{
retrievedCells = ugrid->GetCellAdjacentCells(-1);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
retrievedCells = ugrid->GetCellAdjacentCells(0);
expectedCells = {2, 1, 3};
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells.clear();
retrievedCells = ugrid->GetCellAdjacentCells(4);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells.clear();
retrievedCells = ugrid->GetCellAdjacentCells(5);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells.clear();
retrievedCells.clear();
expectedCells = {1, 5};
retrievedCells = ugrid2d->GetCellAdjacentCells(0);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {0, 2, 3, 4};
retrievedCells = ugrid2d->GetCellAdjacentCells(1);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {1, 3, 4};
retrievedCells = ugrid2d->GetCellAdjacentCells(2);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {2, 1, 4};
retrievedCells = ugrid2d->GetCellAdjacentCells(3);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {1, 2, 3};
retrievedCells = ugrid2d->GetCellAdjacentCells(4);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {0};
retrievedCells = ugrid2d->GetCellAdjacentCells(5);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {1, 5};
retrievedCells = ugrid3d->GetCellAdjacentCells(0);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {0, 2, 5};
retrievedCells = ugrid3d->GetCellAdjacentCells(1);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {1, 4, 3};
retrievedCells = ugrid3d->GetCellAdjacentCells(2);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {4, 2};
retrievedCells = ugrid3d->GetCellAdjacentCells(3);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {2, 3};
retrievedCells = ugrid3d->GetCellAdjacentCells(4);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
expectedCells = {0, 1};
retrievedCells = ugrid3d->GetCellAdjacentCells(5);
TS_ASSERT_EQUALS(expectedCells, retrievedCells);
}
Example - Get a Plan View Polygon
This example shows how to get a plan view polygon. This function first assumes that the shape is prismatic (or that the sides of the shape are vertical) and it gathers the points of the vertical sides. If this method fails, it builds a convex hull of the unique points of the shape. If the second method is employed, concave shapes will result in a convex polygon. The xms::XmUGrid::GetCellPlanViewPolygon function takes a cell index and a vector of Pt3d to contain on return the locations of the points of the plan view polygon and returns whether the operation was succesful. The testing code for this example is XmUGridUnitTests::testGetCellPlanViewPolygon.
{
{{0, 0, 0}, {10, 0, 0}, {10, 10, 0}, {0, 10, 0}},
{{10, 0, 0}, {20, 0, 0}, {20, 10, 0}, {10, 10, 0}},
{{20, 0, 0}, {30, 0, 0}, {20, 10, 0}},
{{30, 0, 0}, {40, 0, 0}, {40, 10, 0}, {40, 20, 0}, {30, 20, 0}, {20, 10, 0}},
{},
{}};
TS_ASSERT(!grid2d->GetCellPlanViewPolygon(-1, viewPolygon));
TS_ASSERT_EQUALS(empty, viewPolygon);
TS_ASSERT(!grid2d->GetCellPlanViewPolygon(grid2d->GetCellCount(), viewPolygon));
TS_ASSERT_EQUALS(empty, viewPolygon);
for (int i = 0; i < grid2d->GetCellCount(); ++i)
{
TS_ASSERT(grid2d->GetCellPlanViewPolygon(i, viewPolygon) != i >= 4);
TS_ASSERT_EQUALS(expectedPolygons[i], viewPolygon);
}
expectedPolygons = {{{0, 0, 0}, {10, 0, 0}, {0, 10, 0}},
{{20, 0, 0}, {20, 10, 0}, {10, 10, 0}, {10, 0, 0}},
{{30, 0, 0}, {30, 10, 0}, {20, 10, 0}, {20, 0, 0}},
{{40, 10, 0}, {40, 20, 0}, {30, 20, 0}, {30, 10, 0}},
{{30, 0, 0}, {40, 0, 0}, {40, 10, 0}, {30, 10, 0}},
{{0, 10, 0}, {10, 10, 0}, {10, 20, 0}, {0, 20, 0}}};
for (int i = 0; i < ugrid3d->GetCellCount(); ++i)
{
TS_ASSERT(ugrid3d->GetCellPlanViewPolygon(i, viewPolygon));
TS_ASSERT_EQUALS(expectedPolygons[i], viewPolygon);
}
expectedPolygons = {{{0.0, 0.0, 0.0}, {20.0, 10.0, 0.0}, {40.0, 0.0, 0.0}, {20.0, 50.0, 0.0}}};
TS_ASSERT(chevronUgrid->GetCellPlanViewPolygon(0, viewPolygon));
TS_ASSERT_EQUALS(expectedPolygons[0], viewPolygon);
}
Example - Get Number of Cell Edges
This example shows how to get the number of cell edges in a cell. The xms::XmUGrid::GetCellEdgeCount function takes a cell index and returns the number of edges as an int. Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::testGetCellEdgeCount
Example Get Cell Edge Adjacent Cells
This example shows how to get all cells which share the specified edge with a cell. The xms::XmUGrid::GetCellEdgeAdjacentCells function takes a cell index and an edge index and returns a vector of cell indices. Functionality is shared between 2d and 3d UGrids. The testing code for this example is shared with other examples, XmUGridUnitTests::testCellEdgeAdjacentCellFunctions
{
VecInt2d expectedCells = {{}, {2}, {1}, {}};
VecInt expected2dCells = {-1, 2, 1, -1};
int adjacentCell;
for (int i(0); i < ugrid->GetCellEdgeCount(0); i++)
{
adjacentCells = ugrid->GetCellEdgeAdjacentCells(0, i);
TS_ASSERT_EQUALS(expectedCells[i], adjacentCells);
adjacentCell = ugrid->GetCell2dEdgeAdjacentCell(0, i);
TS_ASSERT_EQUALS(expected2dCells[i], adjacentCell);
}
VecInt2d expectedCellsFromEdge = {{0, 1}, {0, 2}, {0}};
std::vector<XmEdge> edges = {{1, 4}, {3, 4}, {0, 3}};
adjacentCells = ugrid->GetEdgeAdjacentCells(XmEdge(-1, -1));
for (int i(0); i < edges.size(); i++)
{
adjacentCells = ugrid->GetEdgeAdjacentCells(edges[i]);
TS_ASSERT_EQUALS(expectedCellsFromEdge[i], adjacentCells);
}
}
Example Get Edges Associated with a Cell
This example shows how to get edges associated with a cell. The xms::XmUGrid::GetCellEdges function takes a cell index and returns a vector of XmEdge (each of with contains a pairs of integers which are point indices). Functionality is shared between 2d and 3d UGrids. The testing code for this example is XmUGridUnitTests::testCellEdges
{
std::vector<XmEdge> edges = ugrid->GetCellEdges(0);
std::vector<XmEdge> expectededges = {{0, 3}, {3, 4}, {4, 1}, {1, 0}};
TS_ASSERT_EQUALS(expectededges, edges);
}
Example Get Number of Faces for a Cell
This example shows how to get the number of faces that a cell has. The xms::XmUGrid::GetCell3dFaceCount function takes a cell index and returns the number of cell faces associated with it as an integer. The testing code for this example is XmUGridUnitTests::testCell3dFaceFunctions
{
TS_ASSERT_EQUALS(-1, ugrid2d->GetCell3dFaceCount(-1));
for (int cellIdx = 0; cellIdx < ugrid2d->GetCellCount(); ++cellIdx)
{
TS_ASSERT_EQUALS(0, ugrid2d->GetCell3dFaceCount(cellIdx));
}
TS_ASSERT_EQUALS(-1, ugrid2d->GetCell3dFaceCount(ugrid2d->GetCellCount()));
TS_ASSERT_EQUALS(-1, ugrid3d->GetCell3dFaceCount(-1));
TS_ASSERT_EQUALS(-1, ugrid3d->GetCell3dFaceCount(6));
TS_ASSERT_EQUALS(4, ugrid3d->GetCell3dFaceCount(0));
TS_ASSERT_EQUALS(6, ugrid3d->GetCell3dFaceCount(1));
TS_ASSERT_EQUALS(6, ugrid3d->GetCell3dFaceCount(2));
TS_ASSERT_EQUALS(6, ugrid3d->GetCell3dFaceCount(3));
TS_ASSERT_EQUALS(5, ugrid3d->GetCell3dFaceCount(4));
TS_ASSERT_EQUALS(5, ugrid3d->GetCell3dFaceCount(5));
}
Example Get Cell Face
This example shows how to get a cell face. The xms::XmUGrid::GetCell3dFacePoints function takes a cell index and a face index and returns a vector of point indices which define the face in counter-clockwise direction. The testing code for this example is shared with another example, XmUGridUnitTests::testCell3dFaceFunctions.
Example Get Faces of Cell
This example shows how to get all the faces of a cell. The xms::XmUGrid::GetCell3dFacesPoints function takes a cell index and returns a vector of vectors (one for each face, of point indices which define that face, in counter-clockwise direction as viewed from outside). The testing code for this example is shared with another example, XmUGridUnitTests::testCell3dFaceFunctions.
{
{0, 1, 15},
{1, 5, 15},
{5, 0, 15},
{0, 5, 1},
{1, 16, 21, 6},
{2, 7, 22, 17},
{1, 2, 17, 16},
{6, 21, 22, 7},
{1, 6, 7, 2},
{16, 17, 22, 21},
{2, 17, 22, 7},
{3, 8, 23, 18},
{2, 3, 18, 17},
{7, 22, 23, 8},
{2, 7, 8, 3},
{17, 18, 23, 22},
{9, 8, 13, 14},
{8, 9, 24, 23},
{9, 14, 29, 24},
{14, 13, 28, 29},
{8, 13, 28, 23},
{23, 24, 29, 28},
{3, 4, 18},
{8, 23, 9},
{3, 8, 9, 4},
{4, 9, 23, 18},
{18, 23, 8, 3},
{5, 10, 11, 6},
{5, 6, 20},
{6, 11, 20},
{11, 10, 20},
{10, 5, 20}};
int expectedIdx = 0;
for (int i(0); i < ugrid3d->GetCellCount(); i++)
{
for (int j(0); j < ugrid3d->GetCell3dFaceCount(i); j++, expectedIdx++)
{
cellFace = ugrid3d->GetCell3dFacePoints(i, j);
TS_ASSERT_EQUALS(expectedCellFaces[expectedIdx], cellFace);
}
}
expectedIdx = 0;
for (int i(0); i < ugrid3d->GetCellCount(); i++)
{
cellFaces = ugrid3d->GetCell3dFacesPoints(i);
TS_ASSERT_EQUALS(cellFaces.size(), ugrid3d->GetCell3dFaceCount(i));
for (int j(0); j < ugrid3d->GetCell3dFaceCount(i); j++, expectedIdx++)
{
TS_ASSERT_EQUALS(expectedCellFaces[expectedIdx], cellFaces[j]);
}
}
for (int cellIdx = 0; cellIdx < ugrid3d->GetCellCount(); ++cellIdx)
{
for (int faceIdx = 0; faceIdx < ugrid3d->GetCell3dFaceCount(cellIdx); ++faceIdx)
{
faceOrientations.push_back(ugrid3d->GetCell3dFaceOrientation(cellIdx, faceIdx));
}
}
VecInt expectedFaceOrientations = {
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_TOP, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_BOTTOM,
XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE,
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_TOP,
XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE,
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_TOP,
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_SIDE,
XMU_ORIENTATION_SIDE, XMU_ORIENTATION_TOP,
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_SIDE, XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_TOP,
XMU_ORIENTATION_SIDE,
XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_BOTTOM, XMU_ORIENTATION_TOP, XMU_ORIENTATION_TOP,
XMU_ORIENTATION_SIDE};
TS_ASSERT_EQUALS(expectedFaceOrientations, faceOrientations);
}
Example Get Cell Face Ajacent Cell
This example shows how to get the cell, if any, which shares a face with a specified cell and face. The xms::XmUGrid::GetCell3dFaceAdjacentCell function takes a cell index and a face index and return the cell index of the cell which shares the face, or -1 if there are none. There is an overload which takes a cell index, a face index, an integer which it will set to be the cell index if any of the adjacent cell, and an integer which it will set to be the index of the face of the adjacent cell which it shares with the specified cell if any, and returns true if there is a cell with a shared face. The testing code for this example is XmUGridUnitTests::testGetCell3dFaceAdjacentCell.
{
int rows = 3;
int cols = 2;
int lays = 2;
VecInt expectedNeighbor = {-1, -1, -1, 1, -1, -1, -1, -1, 0, -1, -1, -1};
VecInt expectedFace = {-1, -1, -1, 2, -1, -1, -1, -1, 3, -1, -1, -1};
int expectedIdx = 0;
int neighborCellIdx;
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(-1, 0);
TS_ASSERT_EQUALS(-1, neighborCellIdx);
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(0, -1);
TS_ASSERT_EQUALS(-1, neighborCellIdx);
for (int cellIdx = 0; cellIdx < grid->GetCellCount(); cellIdx++)
{
for (int faceIdx = 0; faceIdx < grid->GetCell3dFaceCount(cellIdx); faceIdx++, expectedIdx++)
{
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(cellIdx, faceIdx);
TS_ASSERT_EQUALS(expectedNeighbor[expectedIdx], neighborCellIdx);
int neighborFaceIdx = -1;
grid->GetCell3dFaceAdjacentCell(cellIdx, faceIdx, neighborCellIdx, neighborFaceIdx);
TS_ASSERT_EQUALS(expectedNeighbor[expectedIdx], neighborCellIdx);
TS_ASSERT_EQUALS(expectedFace[expectedIdx], neighborFaceIdx);
VecInt expectedFacePt = grid->GetCell3dFacePoints(cellIdx, faceIdx);
std::sort(expectedFacePt.begin(), expectedFacePt.end());
VecInt neighborFacePt = grid->GetCell3dFacePoints(neighborCellIdx, neighborFaceIdx);
std::sort(neighborFacePt.begin(), neighborFacePt.end());
if (!expectedFacePt.empty() && !neighborFacePt.empty())
{
TS_ASSERT_EQUALS(expectedFacePt, neighborFacePt);
}
}
}
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(grid->GetCellCount(), 0);
TS_ASSERT_EQUALS(-1, neighborCellIdx);
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(0, grid->GetCellCount());
TS_ASSERT_EQUALS(-1, neighborCellIdx);
rows = 3;
cols = 3;
lays = 3;
expectedNeighbor.clear();
expectedNeighbor = {-1, 1, -1, 2, -1, 4, 0, -1, -1, 3, -1, 5, -1, 3, 0, -1, -1, 6, 2, -1, 1, -1, -1, 7, -1, 5, -1, 6, 0, -1, 4, -1, -1, 7, 1, -1, -1, 7, 4, -1, 2, -1, 6, -1, 5, -1, 3, -1};
expectedFace = {-1, 0, -1, 2, -1, 4, 1, -1, -1, 2, -1, 4, -1, 0, 3, -1, -1, 4, 1, -1, 3, -1, -1, 4, -1, 0, -1, 2, 5, -1, 1, -1, -1, 2, 5, -1, -1, 0, 3, -1, 5, -1, 1, -1, 3, -1, 5, -1};
expectedIdx = 0;
for (int cellIdx = 0; cellIdx < grid->GetCellCount(); ++cellIdx)
{
for (int faceIdx = 0; faceIdx < grid->GetCell3dFaceCount(cellIdx); faceIdx++, expectedIdx++)
{
neighborCellIdx = grid->GetCell3dFaceAdjacentCell(cellIdx, faceIdx);
TS_ASSERT_EQUALS(expectedNeighbor[expectedIdx], neighborCellIdx);
int neighborFaceIdx = -1;
grid->GetCell3dFaceAdjacentCell(cellIdx, faceIdx, neighborCellIdx, neighborFaceIdx);
TS_ASSERT_EQUALS(expectedNeighbor[expectedIdx], neighborCellIdx);
TS_ASSERT_EQUALS(expectedFace[expectedIdx], neighborFaceIdx);
VecInt expectedFacePt = grid->GetCell3dFacePoints(cellIdx, faceIdx);
std::sort(expectedFacePt.begin(), expectedFacePt.end());
VecInt neighborFacePt = grid->GetCell3dFacePoints(neighborCellIdx, neighborFaceIdx);
std::sort(neighborFacePt.begin(), neighborFacePt.end());
if (!expectedFacePt.empty() && !neighborFacePt.empty())
{
TS_ASSERT_EQUALS(expectedFacePt, neighborFacePt);
}
}
}
}