]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Need to add a few more entries to non-local part of sparsity pattern for MeshWorker...
authorkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 20 Jan 2014 09:01:13 +0000 (09:01 +0000)
committerkronbichler <kronbichler@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 20 Jan 2014 09:01:13 +0000 (09:01 +0000)
git-svn-id: https://svn.dealii.org/trunk@32260 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/trilinos_sparse_matrix.h
deal.II/source/dofs/dof_tools_sparsity.cc
deal.II/source/lac/sparsity_tools.cc
deal.II/source/lac/trilinos_sparse_matrix.cc

index 92d4d816f816e98d10a856e6ef0a57f2c48878da..c5ab727b6648044eb76e6755cda9094320f0807b 100644 (file)
@@ -2755,7 +2755,9 @@ namespace TrilinosWrappers
         for (TrilinosWrappers::types::int_type i=0; i<n_columns; ++i)
           std::cout << col_index_ptr[i] << " ";
         std::cout << std::endl << std::endl;
-        std::cout << "Matrix row has the following indices:" << std::endl;
+        std::cout << "Matrix row "
+                  << (row_partitioner().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) == false ? "(nonlocal part)" : "")
+                  << " has the following indices:" << std::endl;
         std::vector<TrilinosWrappers::types::int_type> indices;
         const Epetra_CrsGraph* graph =
           (nonlocal_matrix.get() != 0 &&
index bc26fd13cdd1f10406b03a80c336595025abe02b..6f97a33e9a414fb87730a3d142436d00742b699b 100644 (file)
@@ -596,11 +596,8 @@ namespace DoFTools
                     {
                       // Refinement edges are taken care of by coarser
                       // cells
-
-                      // TODO: in the distributed case, we miss out the
-                      // constraints when the neighbor cell is coarser, but
-                      // only the current cell is owned locally!
-                      if (cell->neighbor_is_coarser(face))
+                      if (cell->neighbor_is_coarser(face) &&
+                          neighbor->subdomain_id() == cell->subdomain_id())
                         continue;
 
                       const unsigned int n_dofs_on_neighbor
@@ -619,14 +616,12 @@ namespace DoFTools
                       // around
                       if (!cell->neighbor(face)->active()
                           ||
-                          (cell->neighbor(face)->subdomain_id() !=
-                           cell->subdomain_id()))
+                          (neighbor->subdomain_id() != cell->subdomain_id()))
                         {
                           constraints.add_entries_local_to_global
                             (dofs_on_other_cell, dofs_on_this_cell,
                              sparsity, keep_constrained_dofs);
-                          if (cell->neighbor(face)->subdomain_id() !=
-                              cell->subdomain_id())
+                          if (neighbor->subdomain_id() != cell->subdomain_id())
                             constraints.add_entries_local_to_global
                             (dofs_on_other_cell, sparsity, keep_constrained_dofs);
                         }
index 57a04430c0718f6ad1571326d82d964925e1fdb3..afa8e0fa3dff3c6bfd79064b2b0148f0ac9a0d2b 100644 (file)
@@ -477,7 +477,7 @@ namespace SparsityTools
           Assert (status.MPI_TAG==124, ExcInternalError());
 
           MPI_Get_count(&status, MPI_BYTE, &len);
-          Assert( len%sizeof(unsigned int)==0, ExcInternalError());
+          Assert( len%sizeof(size_type)==0, ExcInternalError());
 
           recv_buf.resize(len/sizeof(size_type));
 
@@ -500,8 +500,7 @@ namespace SparsityTools
         }
     }
 
-    // complete all sends, so that we can
-    // safely destroy the buffers.
+    // complete all sends, so that we can safely destroy the buffers.
     MPI_Waitall(requests.size(), &requests[0], MPI_STATUSES_IGNORE);
 
   }
@@ -625,8 +624,7 @@ namespace SparsityTools
         }
     }
 
-    // complete all sends, so that we can
-    // safely destroy the buffers.
+    // complete all sends, so that we can safely destroy the buffers.
     MPI_Waitall(requests.size(), &requests[0], MPI_STATUSES_IGNORE);
   }
 
index 55dcb904115f923ed905ad5eb9ea95b821dddff5..88e828f35643f9e7dd88e45f68fcb7485becfcee 100644 (file)
@@ -609,6 +609,9 @@ namespace TrilinosWrappers
                       Utilities::MPI::internal::mpi_type_id(&my_elements),
                       communicator->Comm());
 
+        AssertDimension(std::accumulate(owned_per_proc.begin(),
+          owned_per_proc.end(), size_type()), sparsity_pattern.n_rows());
+
         SparsityTools::distribute_sparsity_pattern
           (const_cast<CompressedSimpleSparsityPattern&>(sparsity_pattern),
            owned_per_proc, communicator->Comm(), sparsity_pattern.row_index_set());
@@ -756,6 +759,7 @@ namespace TrilinosWrappers
 
 
 
+
   void
   SparseMatrix::reinit (const SparsityPattern &sparsity_pattern)
   {
@@ -1804,7 +1808,7 @@ namespace TrilinosWrappers
   template void
   SparseMatrix::reinit (const Epetra_Map &,
                         const Epetra_Map &,
-                        const CompressedSimpleSparsityPattern &,
+                        const CompressedSetSparsityPattern &,
                         const bool);
 
 }

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.