]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Be a bit more carefully when we actually call VecAssembleBegin/End for PETSc vectors.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 3 Aug 2011 15:44:17 +0000 (15:44 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 3 Aug 2011 15:44:17 +0000 (15:44 +0000)
git-svn-id: https://svn.dealii.org/trunk@24003 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/lac/petsc_vector_base.h
deal.II/source/lac/petsc_vector_base.cc

index 08d6c50b3e0bf5d7369deecdda4abc29362cc18e..8f039d9d82259ce3ef97a62301d155479b31db58 100644 (file)
@@ -87,13 +87,6 @@ code using e.g. <code>Tensor<1,dim></code> remains valid.
 <br>
 (Martin Kronbichler, 2011/08/02)
 
-
-<li> Fixed: The function VectorTools::create_right_hand_side now also works
-for objects of type hp::DoFHandler with different finite elements.
-<br>
-(Daniel Gerecht, 2011/07/20)
-
-
 <li> Fixed: deal.II can link with Trilinos but previously it required a
 very specific set of Trilinos sub-libraries; if Trilinos had been compiled
 with a larger set of sub-libraries, linking would sometimes fail. This
@@ -245,6 +238,18 @@ should be fixed now.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> Improved: The PETScWrapper::VectorBase class was rather generous in
+calling the PETSc <code>VecAssembleBegin/End</code> functions that incur
+communication in the %parallel case and are therefore causes of potential
+slowdowns. This has been improved.
+<br>
+(Wolfgang Bangerth, 2011/08/03)
+
+<li> Fixed: The function VectorTools::create_right_hand_side now also works
+for objects of type hp::DoFHandler with different finite elements.
+<br>
+(Daniel Gerecht, 2011/07/20)
+
 <li> Improved: Evaluation of Lagrangian basis functions has been made stable
 by exchanging polynomial evaluation from the standard form
 $a_n x^n+\ldots+a_1 x + a_0$ to a product of linear factors,
index 0f0ca45c2519d62c9526659efb8cb246db6aae23..8bf647553c6be616e6223b4f6c5fab56ec5f6d58 100644 (file)
@@ -906,7 +906,8 @@ namespace PETScWrappers
     const VectorReference &
     VectorReference::operator = (const PetscScalar &value) const
     {
-      if (vector.last_action != VectorBase::LastAction::insert)
+      // if the last action was an addition, we need to flush buffers
+      if (vector.last_action == VectorBase::LastAction::add)
         {
           int ierr;
           ierr = VecAssemblyBegin (vector);
@@ -937,7 +938,8 @@ namespace PETScWrappers
     const VectorReference &
     VectorReference::operator += (const PetscScalar &value) const
     {
-      if (vector.last_action != VectorBase::LastAction::add)
+      // if the last action was a set operation, we need to flush buffers      
+      if (vector.last_action == VectorBase::LastAction::insert)
         {
           int ierr;
           ierr = VecAssemblyBegin (vector);
@@ -945,10 +947,10 @@ namespace PETScWrappers
 
           ierr = VecAssemblyEnd (vector);
           AssertThrow (ierr == 0, ExcPETScError(ierr));
-
-          vector.last_action = VectorBase::LastAction::add;
         }
 
+      vector.last_action = VectorBase::LastAction::add;
+
                                        // we have to do above actions in any
                                        // case to be consistent with the MPI
                                        // communication model (see the
@@ -979,7 +981,8 @@ namespace PETScWrappers
     const VectorReference &
     VectorReference::operator -= (const PetscScalar &value) const
     {
-      if (vector.last_action != VectorBase::LastAction::add)
+      // if the last action was a set operation, we need to flush buffers
+      if (vector.last_action == VectorBase::LastAction::insert)
         {
           int ierr;
           ierr = VecAssemblyBegin (vector);
@@ -987,11 +990,11 @@ namespace PETScWrappers
 
           ierr = VecAssemblyEnd (vector);
           AssertThrow (ierr == 0, ExcPETScError(ierr));
-
-          vector.last_action = VectorBase::LastAction::add;
         }
 
-                                       // we have to do above actions in any
+      vector.last_action = VectorBase::LastAction::add;
+
+                                      // we have to do above actions in any
                                        // case to be consistent with the MPI
                                        // communication model (see the
                                        // comments in the documentation of
@@ -1022,7 +1025,8 @@ namespace PETScWrappers
     const VectorReference &
     VectorReference::operator *= (const PetscScalar &value) const
     {
-      if (vector.last_action != VectorBase::LastAction::insert)
+      // if the last action was an addition, we need to flush buffers
+      if (vector.last_action == VectorBase::LastAction::add)
         {
           int ierr;
           ierr = VecAssemblyBegin (vector);
@@ -1030,10 +1034,10 @@ namespace PETScWrappers
 
           ierr = VecAssemblyEnd (vector);
           AssertThrow (ierr == 0, ExcPETScError(ierr));
-
-          vector.last_action = VectorBase::LastAction::insert;
         }
 
+      vector.last_action = VectorBase::LastAction::insert;
+
                                        // we have to do above actions in any
                                        // case to be consistent with the MPI
                                        // communication model (see the
@@ -1066,7 +1070,8 @@ namespace PETScWrappers
     const VectorReference &
     VectorReference::operator /= (const PetscScalar &value) const
     {
-      if (vector.last_action != VectorBase::LastAction::insert)
+      // if the last action was a set operation, we need to flush buffers
+      if (vector.last_action == VectorBase::LastAction::insert)
         {
           int ierr;
           ierr = VecAssemblyBegin (vector);
@@ -1074,10 +1079,10 @@ namespace PETScWrappers
 
           ierr = VecAssemblyEnd (vector);
           AssertThrow (ierr == 0, ExcPETScError(ierr));
-
-          vector.last_action = VectorBase::LastAction::insert;
         }
 
+      vector.last_action = VectorBase::LastAction::insert;
+
                                        // we have to do above actions in any
                                        // case to be consistent with the MPI
                                        // communication model (see the
index 30b6a33ed97d9272112437ce912a17990af84595..e31f29a01267ed4795cd312a8c10984c54bbae16 100644 (file)
@@ -393,11 +393,22 @@ namespace PETScWrappers
   void
   VectorBase::compress ()
   {
+    // note that one may think that we only need to do something
+    // if in fact the state is anything but LastAction::none. but
+    // that's not true: one frequently gets into situations where
+    // only one processor (or a subset of processors) actually writes
+    // something into a vector, but we still need to call 
+    // VecAssemblyBegin/End on all processors, even those that are
+    // still in LastAction::none state
     int ierr;
     ierr = VecAssemblyBegin(vector);
     AssertThrow (ierr == 0, ExcPETScError(ierr));
     ierr = VecAssemblyEnd(vector);
     AssertThrow (ierr == 0, ExcPETScError(ierr));
+    
+    // reset the last action field to indicate that we're back to
+    // a pristine state
+    last_action = VectorBase::LastAction::none;
   }
 
 
@@ -1110,7 +1121,8 @@ namespace PETScWrappers
                                    const PetscScalar  *values,
                                    const bool          add_values)
   {
-    if (last_action != VectorBase::LastAction::insert)
+    // if the last action was a set operation, we need to flush buffers
+    if (last_action == VectorBase::LastAction::insert)
       {
         int ierr;
         ierr = VecAssemblyBegin (vector);
@@ -1118,6 +1130,10 @@ namespace PETScWrappers
 
         ierr = VecAssemblyEnd (vector);
         AssertThrow (ierr == 0, ExcPETScError(ierr));
+       // reset the last action field to indicate that we're back to
+       // a pristine state
+       last_action = VectorBase::LastAction::none;
       }
 
                                     // VecSetValues complains if we
@@ -1146,6 +1162,8 @@ namespace PETScWrappers
        AssertThrow (ierr == 0, ExcPETScError(ierr));
       }
 
+    // set the mode here, independent of whether we have actually
+    // written elements or whether the list was empty
     last_action = LastAction::insert;
   }
 

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.