]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Attempt a different fix for the problem where local_rows is larger than the total...
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 19 May 2004 14:33:31 +0000 (14:33 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 19 May 2004 14:33:31 +0000 (14:33 +0000)
git-svn-id: https://svn.dealii.org/trunk@9269 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/include/lac/petsc_parallel_sparse_matrix.h
deal.II/lac/source/petsc_parallel_sparse_matrix.cc

index 33c938db9be1513fe3d22df68c03585a39969532..40b1b50f40736666e0fdf940a3c6de95cddd1c93 100644 (file)
@@ -171,8 +171,22 @@ namespace PETScWrappers
                      const std::vector<unsigned int> &row_lengths,
                      const bool                       is_symmetric = false);
 
+                                        /**
+                                         * Exception
+                                         */
+       DeclException2 (ExcLocalRowsTooLarge,
+                       int, int,
+                       << "The number of local rows " << arg1
+                       << " must be larger than the total number of rows " << arg2);
+       
       private:
 
+                                         /**
+                                          * Copy of the communicator object to
+                                          * be used for this parallel vector.
+                                          */
+        MPI_Comm communicator;
+
                                          /**
                                           * Do the actual work for the
                                           * respective reinit() function and the
@@ -194,13 +208,6 @@ namespace PETScWrappers
                         const unsigned int               local_rows,
                         const std::vector<unsigned int> &row_lengths,
                         const bool                       is_symmetric = false);
-
-      private:
-                                         /**
-                                          * Copy of the communicator object to
-                                          * be used for this parallel vector.
-                                          */
-        MPI_Comm communicator;
     };
 
   }
index e71ba0ae962116ae6b9d28d48b3ea14e3675b386..c596fbe0a6ffb1b1a1245fa2232343962abc311e 100644 (file)
@@ -121,6 +121,8 @@ namespace PETScWrappers
                              const unsigned int n_nonzero_per_row,
                              const bool         is_symmetric)
     {
+      Assert (local_rows < m, ExcLocalRowsTooLarge (local_rows, m));
+      
                                        // use the call sequence indicating only
                                        // a maximal number of elements per row
                                        // for all rows globally
@@ -132,7 +134,9 @@ namespace PETScWrappers
                                        // local_rows as they do for local_size
                                        // on parallel vectors
       const int ierr
-        = MatCreateMPIAIJ(communicator, local_rows, local_rows, m, n,
+        = MatCreateMPIAIJ(communicator,
+                         local_rows, std::min(local_rows, n),
+                         m, n,
                           n_nonzero_per_row, 0, 0, 0,
                           &matrix);
       AssertThrow (ierr == 0, ExcPETScError(ierr));
@@ -155,6 +159,8 @@ namespace PETScWrappers
                              const std::vector<unsigned int> &row_lengths,
                              const bool         is_symmetric)
     {
+      Assert (local_rows < m, ExcLocalRowsTooLarge (local_rows, m));
+
       Assert (row_lengths.size() == m,
               ExcDimensionMismatch (row_lengths.size(), m));
     
@@ -176,7 +182,9 @@ namespace PETScWrappers
                                        // local_rows as they do for local_size
                                        // on parallel vectors
       const int ierr
-        = MatCreateMPIAIJ(communicator, local_rows, local_rows, m, n,
+        = MatCreateMPIAIJ(communicator,
+                         local_rows, local_rows,
+                         m, std::min(local_rows, n),
                           0, &int_row_lengths[0], 0, 0,
                           &matrix);
       AssertThrow (ierr == 0, ExcPETScError(ierr));

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.