]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Upgrade PETSc matrices to petsc-dev
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 6 Sep 2011 08:23:59 +0000 (08:23 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 6 Sep 2011 08:23:59 +0000 (08:23 +0000)
git-svn-id: https://svn.dealii.org/trunk@24258 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/petsc_matrix_base.h
deal.II/source/lac/petsc_matrix_base.cc
deal.II/source/lac/petsc_sparse_matrix.cc

index 4e16e9f1f4d0cf65c1727998825aeefd8ad0de2f..5939f5c5ace127ee354a43cbc1eeabca19f81542 100644 (file)
@@ -1108,7 +1108,12 @@ namespace PETScWrappers
                                         * $1000\times32$-bit machine
                                         * precision.
                                         */
-      PetscTruth is_symmetric (const double tolerance = 1.e-12);
+#if DEAL_II_PETSC_VERSION_DEV()
+      PetscBool
+#else
+       PetscTruth 
+#endif
+       is_symmetric (const double tolerance = 1.e-12);
 
 #if DEAL_II_PETSC_VERSION_GTE(2,3,0)
                                         /**
@@ -1120,7 +1125,12 @@ namespace PETScWrappers
                                         * $1000\times32$-bit machine
                                         * precision.
                                         */
-      PetscTruth is_hermitian (const double tolerance = 1.e-12);
+#if DEAL_II_PETSC_VERSION_DEV()
+      PetscBool
+#else
+       PetscTruth 
+#endif
+       is_hermitian (const double tolerance = 1.e-12);
 #endif
 
                                         /*
index d5ca400588b75fab5482e3e2498fb9ea8500c305..ebabbdfd7ab0960fcfb14c846f0136447e0b9fc5 100644 (file)
@@ -90,7 +90,11 @@ namespace PETScWrappers
 
   MatrixBase::~MatrixBase ()
   {
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr = MatDestroy (&matrix);
+#else
     const int ierr = MatDestroy (matrix);
+#endif
     AssertThrow (ierr == 0, ExcPETScError(ierr));
   }
 
@@ -100,7 +104,11 @@ namespace PETScWrappers
   MatrixBase::clear ()
   {
                                      // destroy the matrix...
+#if DEAL_II_PETSC_VERSION_DEV()
+    int ierr = MatDestroy (&matrix);
+#else
     int ierr = MatDestroy (matrix);
+#endif
     AssertThrow (ierr == 0, ExcPETScError(ierr));
                                      // ...and replace it by an empty
                                      // sequential matrix
@@ -146,19 +154,33 @@ namespace PETScWrappers
 #endif
 
     IS index_set;
+#if DEAL_II_PETSC_VERSION_DEV()
+    ISCreateGeneral (get_mpi_communicator(), 1, &petsc_row, PETSC_COPY_VALUES, &index_set);
+#else
     ISCreateGeneral (get_mpi_communicator(), 1, &petsc_row, &index_set);
+#endif
+
 
 #if DEAL_II_PETSC_VERSION_LT(2,3,0)
     const int ierr
       = MatZeroRows(matrix, index_set, &new_diag_value);
+#else
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr
+      = MatZeroRowsIS(matrix, index_set, new_diag_value, PETSC_NULL, PETSC_NULL);
 #else
     const int ierr
       = MatZeroRowsIS(matrix, index_set, new_diag_value);
+#endif
 #endif
 
     AssertThrow (ierr == 0, ExcPETScError(ierr));
 
+#if DEAL_II_PETSC_VERSION_DEV()
+    ISDestroy (&index_set);
+#else
     ISDestroy (index_set);
+#endif
 
     compress ();
   }
@@ -183,20 +205,34 @@ namespace PETScWrappers
                                      // to call them even if #rows is empty,
                                      // since this is a collective operation
     IS index_set;
+#if DEAL_II_PETSC_VERSION_DEV()
+    ISCreateGeneral (get_mpi_communicator(), rows.size(),
+                     &petsc_rows[0], PETSC_COPY_VALUES, &index_set);
+#else
     ISCreateGeneral (get_mpi_communicator(), rows.size(),
                      &petsc_rows[0], &index_set);
+#endif
 
 #if DEAL_II_PETSC_VERSION_LT(2,3,0)
     const int ierr
       = MatZeroRows(matrix, index_set, &new_diag_value);
+#else
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr
+      = MatZeroRowsIS(matrix, index_set, new_diag_value, PETSC_NULL, PETSC_NULL);
 #else
     const int ierr
       = MatZeroRowsIS(matrix, index_set, new_diag_value);
+#endif
 #endif
 
     AssertThrow (ierr == 0, ExcPETScError(ierr));
 
+#if DEAL_II_PETSC_VERSION_DEV()
+    ISDestroy (&index_set);
+#else
     ISDestroy (index_set);
+#endif
 
     compress ();
   }
@@ -570,20 +606,39 @@ namespace PETScWrappers
     AssertThrow (ierr == 0, ExcPETScError(ierr));
   }
 
+#if DEAL_II_PETSC_VERSION_DEV()
+  PetscBool
+#else
   PetscTruth
+#endif
   MatrixBase::is_symmetric (const double tolerance)
   {
-    PetscTruth truth;
+#if DEAL_II_PETSC_VERSION_DEV()
+    PetscBool
+#else
+    PetscTruth 
+#endif
+      truth;
                                        // First flush PETSc caches
     compress ();
     MatIsSymmetric (matrix, tolerance, &truth);
     return truth;
   }
 
+#if DEAL_II_PETSC_VERSION_DEV()
+  PetscBool
+#else
   PetscTruth
+#endif
   MatrixBase::is_hermitian (const double tolerance)
   {
-    PetscTruth truth;
+#if DEAL_II_PETSC_VERSION_DEV()
+    PetscBool
+#else
+    PetscTruth 
+#endif
+      truth;
+
                                     // First flush PETSc caches
     compress ();
 
index 83dbffbd923ccf6dac2190faa1e0dcd546458f34..0fb4b70d5df37bfc93319f8236dd92ad10164b39 100644 (file)
@@ -74,7 +74,7 @@ namespace PETScWrappers
   }
 
 
-
+  
   void
   SparseMatrix::reinit (const unsigned int m,
                         const unsigned int n,
@@ -83,7 +83,11 @@ namespace PETScWrappers
   {
                                      // get rid of old matrix and generate a
                                      // new one
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr = MatDestroy (&matrix);
+#else
     const int ierr = MatDestroy (matrix);
+#endif
     AssertThrow (ierr == 0, ExcPETScError(ierr));
 
     do_reinit (m, n, n_nonzero_per_row, is_symmetric);
@@ -99,7 +103,11 @@ namespace PETScWrappers
   {
                                      // get rid of old matrix and generate a
                                      // new one
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr = MatDestroy (&matrix);
+#else
     const int ierr = MatDestroy (matrix);
+#endif
     AssertThrow (ierr == 0, ExcPETScError(ierr));
 
     do_reinit (m, n, row_lengths, is_symmetric);
@@ -115,7 +123,11 @@ namespace PETScWrappers
   {
                                      // get rid of old matrix and generate a
                                      // new one
+#if DEAL_II_PETSC_VERSION_DEV()
+    const int ierr = MatDestroy (&matrix);
+#else
     const int ierr = MatDestroy (matrix);
+#endif
     AssertThrow (ierr == 0, ExcPETScError(ierr));
 
     do_reinit (sparsity_pattern, preset_nonzero_locations);

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.