[cig-commits] r14388 - in short/3D/PyLith/branches/pylith-swig: libsrc/topology unittests/libtests/topology

brad at geodynamics.org brad at geodynamics.org
Wed Mar 18 21:46:00 PDT 2009


Author: brad
Date: 2009-03-18 21:46:00 -0700 (Wed, 18 Mar 2009)
New Revision: 14388

Modified:
   short/3D/PyLith/branches/pylith-swig/libsrc/topology/Field.cc
   short/3D/PyLith/branches/pylith-swig/unittests/libtests/topology/TestFieldMesh.cc
Log:
Fixed vector scatter for Field.

Modified: short/3D/PyLith/branches/pylith-swig/libsrc/topology/Field.cc
===================================================================
--- short/3D/PyLith/branches/pylith-swig/libsrc/topology/Field.cc	2009-03-19 00:30:11 UTC (rev 14387)
+++ short/3D/PyLith/branches/pylith-swig/libsrc/topology/Field.cc	2009-03-19 04:46:00 UTC (rev 14388)
@@ -472,20 +472,10 @@
   err = VecCreateSeqWithArray(PETSC_COMM_SELF,
 			      _section->sizeWithBC(), _section->restrictSpace(),
 			      &localVec); CHECK_PETSC_ERROR(err);
-  std::cout << "BEFORE LOCALVEC: " << std::endl;
-  VecView(localVec, PETSC_VIEWER_STDOUT_SELF);
-  std::cout << "BEFORE _VECTOR: " << std::endl;
-  VecView(_vector, PETSC_VIEWER_STDOUT_SELF);
-  err = VecScatterBegin(_scatter, _vector, localVec, 
-			INSERT_VALUES, SCATTER_REVERSE); CHECK_PETSC_ERROR(err);
-  err = VecScatterEnd(_scatter, _vector, localVec,
-		      INSERT_VALUES, SCATTER_REVERSE); CHECK_PETSC_ERROR(err);
-
-  std::cout << "AFTER LOCALVEC: " << std::endl;
-  VecView(localVec, PETSC_VIEWER_STDOUT_SELF);
-  std::cout << "AFTER _VECTOR: " << std::endl;
-  VecView(_vector, PETSC_VIEWER_STDOUT_SELF);
-
+  err = VecScatterBegin(_scatter, localVec, _vector,
+			INSERT_VALUES, SCATTER_FORWARD); CHECK_PETSC_ERROR(err);
+  err = VecScatterEnd(_scatter, localVec, _vector,
+		      INSERT_VALUES, SCATTER_FORWARD); CHECK_PETSC_ERROR(err);
   err = VecDestroy(localVec); CHECK_PETSC_ERROR(err);
 } // scatterSectionToVector
 
@@ -508,10 +498,10 @@
   err = VecCreateSeqWithArray(PETSC_COMM_SELF,
 			      _section->sizeWithBC(), _section->restrictSpace(),
 			      &localVec); CHECK_PETSC_ERROR(err);
-  err = VecScatterBegin(_scatter, localVec, _vector,
-			INSERT_VALUES, SCATTER_FORWARD); CHECK_PETSC_ERROR(err);
-  err = VecScatterEnd(_scatter, localVec, _vector,
-		      INSERT_VALUES, SCATTER_FORWARD); CHECK_PETSC_ERROR(err);
+  err = VecScatterBegin(_scatter, _vector, localVec,
+			INSERT_VALUES, SCATTER_REVERSE); CHECK_PETSC_ERROR(err);
+  err = VecScatterEnd(_scatter, _vector, localVec,
+		      INSERT_VALUES, SCATTER_REVERSE); CHECK_PETSC_ERROR(err);
   err = VecDestroy(localVec); CHECK_PETSC_ERROR(err);
 } // scatterVectorToSection
 

Modified: short/3D/PyLith/branches/pylith-swig/unittests/libtests/topology/TestFieldMesh.cc
===================================================================
--- short/3D/PyLith/branches/pylith-swig/unittests/libtests/topology/TestFieldMesh.cc	2009-03-19 00:30:11 UTC (rev 14387)
+++ short/3D/PyLith/branches/pylith-swig/unittests/libtests/topology/TestFieldMesh.cc	2009-03-19 04:46:00 UTC (rev 14388)
@@ -768,11 +768,9 @@
   } // for
 
   field.scatterSectionToVector();
-  CPPUNIT_ASSERT(0 != field._vector);
   CPPUNIT_ASSERT(0 != field._scatter);
   const PetscVec vec = field.vector();
-  std::cout << "VEC: " << std::endl;
-  VecView(vec, PETSC_VIEWER_STDOUT_SELF);
+  CPPUNIT_ASSERT(0 != vec);
   int size = 0;
   VecGetSize(vec, &size);
   double* valuesVec = 0;
@@ -791,6 +789,52 @@
 void
 pylith::topology::TestFieldMesh::testScatterVectorToSection(void)
 { // testScatterVectorToSection
+  const int fiberDim = 3;
+  const double valuesE[] = {
+    1.1, 2.2, 3.3,
+    1.2, 2.3, 3.4,
+    1.3, 2.4, 3.5,
+    1.4, 2.5, 3.6,
+  };
+
+  Mesh mesh;
+  _buildMesh(&mesh);
+  const ALE::Obj<Mesh::SieveMesh>& sieveMesh = mesh.sieveMesh();
+  Field<Mesh> field(mesh);
+  field.newSection(Field<Mesh>::VERTICES_FIELD, fiberDim);
+  field.allocate();
+  const ALE::Obj<Mesh::RealSection>& section = field.section();
+  const ALE::Obj<Mesh::SieveMesh::label_sequence>& vertices = 
+    sieveMesh->depthStratum(0);
+
+  field.createVector();
+  const PetscVec vec = field.vector();
+  CPPUNIT_ASSERT(0 != vec);
+  int size = 0;
+  VecGetSize(vec, &size);
+  double* valuesVec = 0;
+  VecGetArray(vec, &valuesVec);
+
+  const double tolerance = 1.0e-06;
+  const int sizeE = vertices->size() * fiberDim;
+  CPPUNIT_ASSERT_EQUAL(sizeE, size);
+  for (int i=0; i < sizeE; ++i)
+    valuesVec[i] = valuesE[i];
+  VecRestoreArray(vec, &valuesVec);
+
+  field.scatterVectorToSection();
+  CPPUNIT_ASSERT(0 != field._scatter);
+
+  double_array values(fiberDim);
+  int i = 0;
+  for (Mesh::SieveMesh::label_sequence::iterator v_iter=vertices->begin();
+       v_iter != vertices->end();
+       ++v_iter) {
+    section->restrictPoint(*v_iter, &values[0], fiberDim);
+    for (int iDim=0; iDim < fiberDim; ++iDim)
+      CPPUNIT_ASSERT_DOUBLES_EQUAL(valuesE[i++], values[iDim], tolerance);
+  } // for
+
 } // testScatterVectorToSection
 
 // ----------------------------------------------------------------------



More information about the CIG-COMMITS mailing list