[cig-commits] r6974 - in short/3D/PyLith/trunk: . modulesrc/utils pylith/problems

brad at geodynamics.org brad at geodynamics.org
Fri May 25 17:04:04 PDT 2007


Author: brad
Date: 2007-05-25 17:04:03 -0700 (Fri, 25 May 2007)
New Revision: 6974

Modified:
   short/3D/PyLith/trunk/TODO
   short/3D/PyLith/trunk/modulesrc/utils/petsc.pyxe.src
   short/3D/PyLith/trunk/pylith/problems/Explicit.py
   short/3D/PyLith/trunk/pylith/problems/Implicit.py
Log:
Added bindings for PETSc MatZeroEntries().

Modified: short/3D/PyLith/trunk/TODO
===================================================================
--- short/3D/PyLith/trunk/TODO	2007-05-25 23:31:49 UTC (rev 6973)
+++ short/3D/PyLith/trunk/TODO	2007-05-26 00:04:03 UTC (rev 6974)
@@ -2,9 +2,6 @@
 MAIN PRIORITIES (Brad)
 ======================================================================
 
-0. Create binding for zeroing PETSc matrix (use in Explicit.py and
-     Implicit.py before integrating Jacobian)
-
 1. Unit tests for Integrator stuff.
 
   a. Add in use of useElasticBehavior()

Modified: short/3D/PyLith/trunk/modulesrc/utils/petsc.pyxe.src
===================================================================
--- short/3D/PyLith/trunk/modulesrc/utils/petsc.pyxe.src	2007-05-25 23:31:49 UTC (rev 6973)
+++ short/3D/PyLith/trunk/modulesrc/utils/petsc.pyxe.src	2007-05-26 00:04:03 UTC (rev 6974)
@@ -103,4 +103,23 @@
   return
 
 
+def mat_setzero(mat):
+  """
+  Zero out entries in matrix (retain structure).
+  """
+  # create shim for 'MatZeroEntries'
+  #embed{ int Mat_setzero(void* matVptr)
+  Mat* mat = (Mat*) matVptr;
+  PetscErrorCode err = 0;
+  err = MatZeroEntries(*mat);
+  return 0;
+  #}embed
+  cdef void* matVptr
+  matVptr = PyCObject_AsVoidPtr(mat)
+  err = Mat_setzero(matVptr)
+  if err:
+    raise RuntimError("Error zeroing matrix values.")
+  return
+
+
 # End of file 

Modified: short/3D/PyLith/trunk/pylith/problems/Explicit.py
===================================================================
--- short/3D/PyLith/trunk/pylith/problems/Explicit.py	2007-05-25 23:31:49 UTC (rev 6973)
+++ short/3D/PyLith/trunk/pylith/problems/Explicit.py	2007-05-26 00:04:03 UTC (rev 6974)
@@ -81,7 +81,7 @@
 
     self._info.log("Forming Jacobian of operator.")
     import pylith.utils.petsc as petsc
-    #petsc.zeroMatrix(self.jacobian)
+    petsc.mat_setzero(self.jacobian)
     for integrator in self.integrators:
       integrator.timeStep(dt)
       integrator.integrateJacobian(self.jacobian, self.fields.cppHandle)
@@ -116,7 +116,7 @@
     if needNewJacobian:
       self._info.log("Reforming Jacobian of operator.")
       import pylith.utils.petsc as petsc
-      petsc.zeroMatrix(self.jacobian)
+      petsc.mat_setzero(self.jacobian)
       for integrator in self.integrators:
         integrator.timeStep(dt)
         integrator.integrateJacobian(self.jacobian, self.fields.cppHandle)

Modified: short/3D/PyLith/trunk/pylith/problems/Implicit.py
===================================================================
--- short/3D/PyLith/trunk/pylith/problems/Implicit.py	2007-05-25 23:31:49 UTC (rev 6973)
+++ short/3D/PyLith/trunk/pylith/problems/Implicit.py	2007-05-26 00:04:03 UTC (rev 6974)
@@ -136,7 +136,7 @@
     if needNewJacobian:
       self._info.log("Reforming Jacobian of operator.")
       import pylith.utils.petsc as petsc
-      petsc.zeroMatrix(self.jacobian)
+      petsc.mat_setzero(self.jacobian)
       for integrator in self.integrators:
         integrator.timeStep(dt)
         integrator.integrateJacobian(self.jacobian, self.fields.cppHandle)
@@ -201,10 +201,14 @@
 
     self._info.log("Setting constraints.")
     dispT = self.fields.getReal("dispT")
+    import pylith.topology.topology as bindings
+    bindings.zeroRealSection(dispT)
     for constraint in self.constraints:
       constraint.setField(t, dispT)
 
     self._info.log("Integrating Jacobian and residual of operator.")
+    import pylith.utils.petsc as petsc
+    petsc.mat_setzero(self.jacobian)
     for integrator in self.integrators:
       integrator.timeStep(dt)
       integrator.integrateJacobian(self.jacobian, self.fields.cppHandle)



More information about the cig-commits mailing list