[cig-commits] r5119 - in short/3D/PyLith/trunk: libsrc/feassemble unittests/libtests/feassemble

knepley at geodynamics.org knepley at geodynamics.org
Mon Oct 30 14:48:56 PST 2006


Author: knepley
Date: 2006-10-30 14:48:56 -0800 (Mon, 30 Oct 2006)
New Revision: 5119

Modified:
   short/3D/PyLith/trunk/libsrc/feassemble/Integrator.hh
   short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.cc
   short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.hh
   short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.cc
   short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.hh
   short/3D/PyLith/trunk/unittests/libtests/feassemble/TestIntegrator.cc
Log:
Fixed compliation/link problems


Modified: short/3D/PyLith/trunk/libsrc/feassemble/Integrator.hh
===================================================================
--- short/3D/PyLith/trunk/libsrc/feassemble/Integrator.hh	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/libsrc/feassemble/Integrator.hh	2006-10-30 22:48:56 UTC (rev 5119)
@@ -75,7 +75,6 @@
    */
   virtual 
   void integrate(PetscMat* mat,
-		 ALE::Obj<ALE::Mesh>& mesh,
 		 const ALE::Obj<real_section_type>& fieldIn,
 		 const ALE::Obj<real_section_type>& coordinates) = 0;
 

Modified: short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.cc
===================================================================
--- short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.cc	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.cc	2006-10-30 22:48:56 UTC (rev 5119)
@@ -50,7 +50,6 @@
 // Compute matrix associated with operator.
 void
 pylith::feassemble::IntegratorElasticity3D::integrate(PetscMat* mat,
-		    ALE::Obj<ALE::Mesh>& mesh,
    		    const ALE::Obj<ALE::Mesh::real_section_type>& fieldIn,
 		    const ALE::Obj<ALE::Mesh::real_section_type>& coordinates)
 { // integrate

Modified: short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.hh
===================================================================
--- short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.hh	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/libsrc/feassemble/IntegratorElasticity3D.hh	2006-10-30 22:48:56 UTC (rev 5119)
@@ -60,7 +60,6 @@
    * @param coordinates Field of cell vertex coordinates
    */
   void integrate(PetscMat* mat,
-		 ALE::Obj<ALE::Mesh>& mesh,
 		 const ALE::Obj<ALE::Mesh::real_section_type>& fieldIn,
 		 const ALE::Obj<ALE::Mesh::real_section_type>& coordinates);
 

Modified: short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.cc
===================================================================
--- short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.cc	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.cc	2006-10-30 22:48:56 UTC (rev 5119)
@@ -115,7 +115,6 @@
 void
 pylith::feassemble::IntegratorInertia::integrate(
 			     PetscMat* mat,
-			     ALE::Obj<ALE::Mesh>& mesh,
 			     const ALE::Obj<real_section_type>& fieldIn,
 			     const ALE::Obj<real_section_type>& coordinates)
 { // integrate
@@ -142,7 +141,7 @@
   err = MatCreate(topology->comm(), mat);
   err = MatSetSizes(*mat, localSize, localSize, globalSize, globalSize);
   err = MatSetFromOptions(*mat);
-  err = preallocateMatrix(mesh, fieldIn->getAtlas(), globalOrder, *mat);
+  err = preallocateMatrix(topology, fieldIn->getAtlas(), globalOrder, *mat);
 
   // Allocate matrix for cell values (if necessary)
   _initCellMatrix();

Modified: short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.hh
===================================================================
--- short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.hh	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/libsrc/feassemble/IntegratorInertia.hh	2006-10-30 22:48:56 UTC (rev 5119)
@@ -75,7 +75,6 @@
    * @param coordinates Field of cell vertex coordinates
    */
   void integrate(PetscMat* mat,
-		 ALE::Obj<ALE::Mesh>& mesh,
                  const ALE::Obj<real_section_type>& fieldIn,
                  const ALE::Obj<real_section_type>& coordinates);
 

Modified: short/3D/PyLith/trunk/unittests/libtests/feassemble/TestIntegrator.cc
===================================================================
--- short/3D/PyLith/trunk/unittests/libtests/feassemble/TestIntegrator.cc	2006-10-30 21:35:42 UTC (rev 5118)
+++ short/3D/PyLith/trunk/unittests/libtests/feassemble/TestIntegrator.cc	2006-10-30 22:48:56 UTC (rev 5119)
@@ -147,6 +147,7 @@
 					  const IntegratorData& data) const
 { // _testIntegrate
   typedef ALE::Mesh::real_section_type real_section_type;
+  typedef ALE::Mesh::topology_type topology_type;
 
   ALE::Obj<ALE::Mesh> mesh = _TestIntegrator::_setupMesh(data);
   const ALE::Mesh::int_section_type::patch_type patch = 0;
@@ -160,14 +161,23 @@
   fieldIn->setName("fieldIn");
   fieldIn->setFiberDimensionByDepth(patch, 0, fiberDim);
   fieldIn->allocate();
+  int iVertex = 0;
+  const ALE::Obj<topology_type::label_sequence>& vertices = 
+    mesh->getTopology()->depthStratum(patch, 0);
+  const topology_type::label_sequence::iterator verticesEnd =
+    vertices->end();
+  for (topology_type::label_sequence::iterator vIter=vertices->begin();
+       vIter != verticesEnd;
+       ++vIter, ++iVertex)
+    fieldIn->update(patch, *vIter, &data.fieldIn[iVertex*fiberDim]);
 
   // Integrate
   PetscMat mat;
   const ALE::Obj<real_section_type>& coordinates = 
     mesh->getRealSection("coordinates");
-  integrator->integrate(&mat, mesh, fieldIn, coordinates);
+  integrator->integrate(&mat, fieldIn, coordinates);
 
-  // Crate dense matrix
+  // Create dense matrix
   // :TODO: ADD STUFF HERE
 
   // Get values associated with dense matrix



More information about the cig-commits mailing list