[cig-commits] r21439 - in short/3D/PyLith/trunk/unittests/libtests/bc: . data

knepley at geodynamics.org knepley at geodynamics.org
Sun Mar 3 13:32:02 PST 2013


Author: knepley
Date: 2013-03-03 13:32:02 -0800 (Sun, 03 Mar 2013)
New Revision: 21439

Modified:
   short/3D/PyLith/trunk/unittests/libtests/bc/TestBoundaryMesh.cc
   short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataHex8.cc
   short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataQuad4.cc
   short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTet4.cc
   short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTri3.cc
   short/3D/PyLith/trunk/unittests/libtests/bc/data/hex8.mesh
   short/3D/PyLith/trunk/unittests/libtests/bc/data/quad4.mesh
   short/3D/PyLith/trunk/unittests/libtests/bc/data/tri3.mesh
Log:
Fixed boundary mesh tests for DMPlex

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/TestBoundaryMesh.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/TestBoundaryMesh.cc	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/TestBoundaryMesh.cc	2013-03-03 21:32:02 UTC (rev 21439)
@@ -57,6 +57,7 @@
   CPPUNIT_ASSERT(0 != _data);
 
   topology::Mesh mesh;
+  PetscErrorCode err;
 
   meshio::MeshIOAscii iohandler;
   iohandler.filename(_data->filename);
@@ -72,45 +73,42 @@
 
   // Create submesh
   topology::SubMesh submesh(mesh, _data->bcLabel);
-  CPPUNIT_ASSERT(submesh.dmMesh());
+  DM                dmMesh = submesh.dmMesh();
+  CPPUNIT_ASSERT(dmMesh);
 
   // Check vertices
-  const ALE::Obj<SieveSubMesh::label_sequence>& vertices = 
-    submesh.sieveMesh()->depthStratum(0);
-  const SieveSubMesh::label_sequence::iterator verticesEnd = vertices->end();
-  CPPUNIT_ASSERT_EQUAL(_data->numVerticesNoFault, int(vertices->size()));
+  IS              subpointIS;
+  const PetscInt *subpointMap;
+  PetscInt        vStart, vEnd;
 
-  int ipt = 0;
-  for (SieveSubMesh::label_sequence::iterator v_iter=vertices->begin();
-       v_iter != verticesEnd;
-       ++v_iter, ++ipt)
-    CPPUNIT_ASSERT_EQUAL(_data->verticesNoFault[ipt], *v_iter);
+  err = DMPlexGetDepthStratum(dmMesh, 0, &vStart, &vEnd);CHECK_PETSC_ERROR(err);
+  CPPUNIT_ASSERT_EQUAL(_data->numVerticesNoFault, vEnd-vStart);
+  err = DMPlexCreateSubpointIS(dmMesh, &subpointIS);CHECK_PETSC_ERROR(err);
+  err = ISGetIndices(subpointIS, &subpointMap);CHECK_PETSC_ERROR(err);
+  for (PetscInt v = vStart; v < vEnd; ++v)
+    CPPUNIT_ASSERT_EQUAL(_data->verticesNoFault[v-vStart], subpointMap[v]);
 
   // Check cells
-  const ALE::Obj<SieveSubMesh::label_sequence>& cells = 
-    submesh.sieveMesh()->heightStratum(1);
-  const SieveSubMesh::label_sequence::iterator cellsEnd = cells->end();
-  const ALE::Obj<SieveSubMesh::sieve_type>& sieve = 
-    submesh.sieveMesh()->getSieve();
-  assert(!sieve.isNull());
-  CPPUNIT_ASSERT_EQUAL(_data->numCells, int(cells->size()));
+  PetscInt cStart, cEnd;
 
-  ALE::ISieveVisitor::NConeRetriever<SieveSubMesh::sieve_type> ncV(*sieve, (int) pow(sieve->getMaxConeSize(), submesh.sieveMesh()->depth()));
+  err = DMPlexGetHeightStratum(dmMesh, 1, &cStart, &cEnd);CHECK_PETSC_ERROR(err);
+  CPPUNIT_ASSERT_EQUAL(_data->numCells, cEnd-cStart);
 
-  int icell = 0;
-  int index = 0;
-  for (SieveSubMesh::label_sequence::iterator c_iter=cells->begin();
-       c_iter != cellsEnd;
-       ++c_iter, ++icell) {
-    ALE::ISieveTraversal<SieveSubMesh::sieve_type>::orientedClosure(*sieve, *c_iter, ncV);
-    const int coneSize = ncV.getSize();
-    const SieveSubMesh::point_type *cone = ncV.getPoints();
-    CPPUNIT_ASSERT_EQUAL(_data->numCorners, coneSize);
+  for (PetscInt c = cStart, index = 0; c < cEnd; ++c) {
+    PetscInt *closure = NULL;
+    PetscInt  closureSize, numVertices = 0;
 
-    for(int v = 0; v < coneSize; ++v, ++index)
-      CPPUNIT_ASSERT_EQUAL(_data->cellsNoFault[index], cone[v]);
-    ncV.clear();
+    err = DMPlexGetTransitiveClosure(dmMesh, c, PETSC_TRUE, &closureSize, &closure);CHECK_PETSC_ERROR(err);
+    for (PetscInt p = 0; p < closureSize*2; p += 2) {
+      if ((closure[p] >= vStart) && (closure[p] < vEnd)) closure[numVertices++] = closure[p];
+    }
+    CPPUNIT_ASSERT_EQUAL(_data->numCorners, numVertices);
+    for (PetscInt v = 0; v < numVertices; ++v, ++index)
+      CPPUNIT_ASSERT_EQUAL(_data->cellsNoFault[index], subpointMap[closure[v]]);
+    err = DMPlexRestoreTransitiveClosure(dmMesh, c, PETSC_TRUE, &closureSize, &closure);CHECK_PETSC_ERROR(err);
   } // for
+  err = ISRestoreIndices(subpointIS, &subpointMap);CHECK_PETSC_ERROR(err);
+  err = ISDestroy(&subpointIS);CHECK_PETSC_ERROR(err);
 } // testSubmesh
 
 // ----------------------------------------------------------------------
@@ -119,6 +117,7 @@
 pylith::bc::TestBoundaryMesh::testSubmeshFault(void)
 { // testSubmeshFault
   CPPUNIT_ASSERT(0 != _data);
+  PetscErrorCode err;
 
   topology::Mesh mesh;
   meshio::MeshIOAscii iohandler;
@@ -135,54 +134,52 @@
 
   // Adjust topology
   faults::FaultCohesiveKin fault;
-  int firstFaultVertex    = 0;
-  int firstLagrangeVertex = mesh.sieveMesh()->getIntSection(_data->faultLabel)->size();
-  int firstFaultCell      = mesh.sieveMesh()->getIntSection(_data->faultLabel)->size();
+  PetscInt firstFaultVertex = 0;
+  PetscInt firstLagrangeVertex, firstFaultCell;
+
+  err = DMPlexGetStratumSize(mesh.dmMesh(), _data->faultLabel, 1, &firstLagrangeVertex);CHECK_PETSC_ERROR(err);
+  firstFaultCell = firstLagrangeVertex;
   fault.label(_data->faultLabel);
   fault.id(_data->faultId);
   fault.adjustTopology(&mesh, &firstFaultVertex, &firstLagrangeVertex, &firstFaultCell, _flipFault);
 
   // Create submesh
   topology::SubMesh submesh(mesh, _data->bcLabel);
-  CPPUNIT_ASSERT(!submesh.sieveMesh().isNull());
+  DM                dmMesh = submesh.dmMesh();
+  CPPUNIT_ASSERT(dmMesh);
 
   // Check vertices
-  const ALE::Obj<SieveSubMesh::label_sequence>& vertices = 
-    submesh.sieveMesh()->depthStratum(0);
-  const SieveSubMesh::label_sequence::iterator verticesEnd = vertices->end();
-  CPPUNIT_ASSERT_EQUAL(_data->numVerticesFault, int(vertices->size()));
+  IS              subpointIS;
+  const PetscInt *subpointMap;
+  PetscInt        vStart, vEnd;
 
-  int ipt = 0;
-  for (SieveSubMesh::label_sequence::iterator v_iter=vertices->begin();
-       v_iter != verticesEnd;
-       ++v_iter, ++ipt)
-    CPPUNIT_ASSERT_EQUAL(_data->verticesFault[ipt], *v_iter);
-
+  err = DMPlexGetDepthStratum(dmMesh, 0, &vStart, &vEnd);CHECK_PETSC_ERROR(err);
+  CPPUNIT_ASSERT_EQUAL(_data->numVerticesFault, vEnd-vStart);
+  err = DMPlexCreateSubpointIS(dmMesh, &subpointIS);CHECK_PETSC_ERROR(err);
+  err = ISGetIndices(subpointIS, &subpointMap);CHECK_PETSC_ERROR(err);
+  for (PetscInt v = vStart; v < vEnd; ++v)
+    CPPUNIT_ASSERT_EQUAL(_data->verticesFault[v-vStart], subpointMap[v]);
   // Check cells
-  const ALE::Obj<SieveSubMesh::label_sequence>& cells = 
-    submesh.sieveMesh()->heightStratum(1);
-  const SieveSubMesh::label_sequence::iterator cellsEnd = cells->end();
-  const ALE::Obj<SieveSubMesh::sieve_type>& sieve = 
-    submesh.sieveMesh()->getSieve();
-  assert(!sieve.isNull());
-  CPPUNIT_ASSERT_EQUAL(_data->numCells, int(cells->size()));
+  PetscInt cStart, cEnd;
 
-  ALE::ISieveVisitor::NConeRetriever<SieveSubMesh::sieve_type> ncV(*sieve, (int) pow(sieve->getMaxConeSize(), submesh.sieveMesh()->depth()));
+  err = DMPlexGetHeightStratum(dmMesh, 1, &cStart, &cEnd);CHECK_PETSC_ERROR(err);
+  CPPUNIT_ASSERT_EQUAL(_data->numCells, cEnd-cStart);
 
-  int icell = 0;
-  int index = 0;
-  for (SieveSubMesh::label_sequence::iterator c_iter=cells->begin();
-       c_iter != cellsEnd;
-       ++c_iter, ++icell) {
-    ALE::ISieveTraversal<SieveSubMesh::sieve_type>::orientedClosure(*sieve, *c_iter, ncV);
-    const int coneSize = ncV.getSize();
-    const SieveSubMesh::point_type *cone = ncV.getPoints();
-    CPPUNIT_ASSERT_EQUAL(_data->numCorners, coneSize);
+  for (PetscInt c = cStart, index = 0; c < cEnd; ++c) {
+    PetscInt *closure = NULL;
+    PetscInt  closureSize, numVertices = 0;
 
-    for(int v = 0; v < coneSize; ++v, ++index)
-      CPPUNIT_ASSERT_EQUAL(_data->cellsFault[index], cone[v]);
-    ncV.clear();
+    err = DMPlexGetTransitiveClosure(dmMesh, c, PETSC_TRUE, &closureSize, &closure);CHECK_PETSC_ERROR(err);
+    for (PetscInt p = 0; p < closureSize*2; p += 2) {
+      if ((closure[p] >= vStart) && (closure[p] < vEnd)) closure[numVertices++] = closure[p];
+    }
+    CPPUNIT_ASSERT_EQUAL(_data->numCorners, numVertices);
+    for (PetscInt v = 0; v < numVertices; ++v, ++index)
+      CPPUNIT_ASSERT_EQUAL(_data->cellsFault[index], subpointMap[closure[v]]);
+    err = DMPlexRestoreTransitiveClosure(dmMesh, c, PETSC_TRUE, &closureSize, &closure);CHECK_PETSC_ERROR(err);
   } // for
+  err = ISRestoreIndices(subpointIS, &subpointMap);CHECK_PETSC_ERROR(err);
+  err = ISDestroy(&subpointIS);CHECK_PETSC_ERROR(err);
 } // testSubmeshFault
 
 

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataHex8.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataHex8.cc	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataHex8.cc	2013-03-03 21:32:02 UTC (rev 21439)
@@ -39,11 +39,11 @@
 
 const int pylith::bc::BoundaryMeshDataHex8::_numVerticesFault = 8;
 const int pylith::bc::BoundaryMeshDataHex8::_verticesFault[] = {
-  2, 4, 6, 8, 10, 12, 14, 16
+  3, 5, 7, 9, 11, 13, 15, 17
 };
 const int pylith::bc::BoundaryMeshDataHex8::_cellsFault[] = {
-  2, 14, 16, 8,
-  4, 6, 12, 10,
+  3, 15, 17, 9,
+  5, 7, 13, 11,
 };
 
 pylith::bc::BoundaryMeshDataHex8::BoundaryMeshDataHex8(void)

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataQuad4.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataQuad4.cc	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataQuad4.cc	2013-03-03 21:32:02 UTC (rev 21439)
@@ -39,11 +39,11 @@
 
 const int pylith::bc::BoundaryMeshDataQuad4::_numVerticesFault = 4;
 const int pylith::bc::BoundaryMeshDataQuad4::_verticesFault[] = {
-  2, 4, 6, 8,
+  3, 5, 7, 9,
 };
 const int pylith::bc::BoundaryMeshDataQuad4::_cellsFault[] = {
-  2, 8,
-  4, 6,
+  3, 9,
+  5, 7,
 };
 
 pylith::bc::BoundaryMeshDataQuad4::BoundaryMeshDataQuad4(void)

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTet4.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTet4.cc	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTet4.cc	2013-03-03 21:32:02 UTC (rev 21439)
@@ -39,11 +39,11 @@
 
 const int pylith::bc::BoundaryMeshDataTet4::_numVerticesFault = 6;
 const int pylith::bc::BoundaryMeshDataTet4::_verticesFault[] = {
-  2, 4, 5, 6, 7, 9
+  3, 5, 6, 7, 8, 10
 };
 const int pylith::bc::BoundaryMeshDataTet4::_cellsFault[] = {
-  5, 4, 2,
-  7, 6, 9,
+  6, 5, 3,
+  8, 7, 10,
 };
 
 pylith::bc::BoundaryMeshDataTet4::BoundaryMeshDataTet4(void)

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTri3.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTri3.cc	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/BoundaryMeshDataTri3.cc	2013-03-03 21:32:02 UTC (rev 21439)
@@ -38,10 +38,10 @@
 
 const int pylith::bc::BoundaryMeshDataTri3::_numVerticesFault = 2;
 const int pylith::bc::BoundaryMeshDataTri3::_verticesFault[] = {
-  3, 5
+  4, 6
 };
 const int pylith::bc::BoundaryMeshDataTri3::_cellsFault[] = {
-  3, 5
+  4, 6
 };
 
 pylith::bc::BoundaryMeshDataTri3::BoundaryMeshDataTri3(void)

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/hex8.mesh
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/hex8.mesh	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/hex8.mesh	2013-03-03 21:32:02 UTC (rev 21439)
@@ -27,6 +27,24 @@
 // |/      |/    0
 // 3 ----- 2
 //
+//
+// Sieve mesh with fault
+//
+//            14 -----13
+//           / |     / |
+//          /  |    /  |
+//        12 -----11   |
+//       / |  8  /-|-- 7
+//      /  | /  /  | /
+//    18 -----17   |/  1
+//   / |   6 /-|-- 5
+//  /  | /  /  | /
+//10 ----- 9   |/    2
+// |  16 --|--15
+// | /     | /
+// |/      |/    0
+// 4 ----- 3
+//
 mesh = {
   dimension = 3
   use-index-zero = true

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/quad4.mesh
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/quad4.mesh	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/quad4.mesh	2013-03-03 21:32:02 UTC (rev 21439)
@@ -16,11 +16,11 @@
 //
 // Sieve mesh w/fault
 //
-//  3 ----- 5-11-9 ----- 7
+//  4 -----10-12-6 ----- 8
 //  |       |    |       |
-//  |   0   |    |   1   |
+//  |   0   |  2 |   1   |
 //  |       |    |       |
-//  2 ----- 4-10-8 ----- 6
+//  3 ----- 9-11-5 ----- 7
 //
 mesh = {
   dimension = 2

Modified: short/3D/PyLith/trunk/unittests/libtests/bc/data/tri3.mesh
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/bc/data/tri3.mesh	2013-03-02 08:46:07 UTC (rev 21438)
+++ short/3D/PyLith/trunk/unittests/libtests/bc/data/tri3.mesh	2013-03-03 21:32:02 UTC (rev 21439)
@@ -16,6 +16,15 @@
 //    \ | /
 //      3
 //
+// Sieve mesh with fault
+//      8---5
+//    / |   | \
+//   /  |   |  \
+//  3 0 | 2 | 1 6
+//   \  |   |  /
+//    \ |   | /
+//      7---4
+//
 mesh = {
   dimension = 2
   use-index-zero = true



More information about the CIG-COMMITS mailing list