]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Patch by Fahad Alrashed: Overload PETScWrappers::MPI::Vector::print to make sure...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 2 Jan 2013 16:59:56 +0000 (16:59 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 2 Jan 2013 16:59:56 +0000 (16:59 +0000)
git-svn-id: https://svn.dealii.org/trunk@27897 0785d39b-7218-0410-832d-ea1e28bc413d

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

index 486e97a696ccb59d37d1e24433a3f5fb16c557b2..2c4dd85becb6b947ea59d8f1a83ac53286935225 100644 (file)
@@ -127,6 +127,12 @@ DoFHandler, in particular removal of specializations.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> New: The PETScWrappers::MPI::Vector::print function overloads the
+function of same name in the base class to ensure that the output
+generated by a parallel vector makes sense.
+<br>
+(Fahad Alrashed, 2013/1/2)
+
 <li> New: The PETScWrappers::VectorBase class now has a function
 PETScWrappers::VectorBase::write_ascii() that allows writing the
 vector's data to the default output stream.
index 00246610f39e67c32728e6296c97d5f613a7a99d..bea148060f811a1f7664b0bc89c418b16c3a5091 100644 (file)
@@ -1,7 +1,7 @@
 //---------------------------------------------------------------------------
 //    $Id$
 //
-//    Copyright (C) 2004, 2005, 2006, 2007, 2009, 2010, 2012 by the deal.II authors
+//    Copyright (C) 2004, 2005, 2006, 2007, 2009, 2010, 2012, 2013 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -391,6 +391,30 @@ namespace PETScWrappers
        */
       const MPI_Comm &get_mpi_communicator () const;
 
+      /**
+       * Print to a
+       * stream. @p precision denotes
+       * the desired precision with
+       * which values shall be printed,
+       * @p scientific whether
+       * scientific notation shall be
+       * used. If @p across is
+       * @p true then the vector is
+       * printed in a line, while if
+       * @p false then the elements
+       * are printed on a separate line
+       * each.
+       *
+       * @note This function overloads the one in the
+       * base class to ensure that the right thing happens
+       * for parallel vectors that are distributed across
+       * processors.
+       */
+      void print (std::ostream       &out,
+                  const unsigned int  precision  = 3,
+                  const bool          scientific = true,
+                  const bool          across     = true) const;
+
     protected:
       /**
        * Create a vector of length
index 225c1ebda18ef4016d881fdc45616d564113640f..612be271935b63ac4a440a83b9a9652b11087512 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 2004, 2006, 2008, 2009, 2010, 2011, 2012 by the deal.II authors
+//    Copyright (C) 2004, 2006, 2008, 2009, 2010, 2011, 2012, 2013 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -286,6 +286,76 @@ namespace PETScWrappers
 
 
 
+    void
+    Vector::print (std::ostream      &out,
+                   const unsigned int precision,
+                   const bool         scientific,
+                   const bool         across) const
+    {
+      AssertThrow (out, ExcIO());
+
+      // get a representation of the vector and
+      // loop over all the elements
+      PetscScalar *val;
+      PetscInt    nlocal, istart, iend;
+
+      int ierr = VecGetArray (vector, &val);
+
+      AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+      ierr = VecGetLocalSize (vector, &nlocal);
+
+      AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+      ierr = VecGetOwnershipRange (vector, &istart, &iend);
+
+      AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+      // save the state of out stream
+      std::ios::fmtflags old_flags = out.flags();
+      unsigned int old_precision = out.precision (precision);
+
+      out.precision (precision);
+      if (scientific)
+        out.setf (std::ios::scientific, std::ios::floatfield);
+      else
+        out.setf (std::ios::fixed, std::ios::floatfield);
+
+      for ( unsigned int i = 0;
+           i < Utilities::MPI::n_mpi_processes(communicator);
+           i++)
+       {
+         // This is slow, but most likely only used to debug.
+         MPI_Barrier(communicator);
+         if (i == Utilities::MPI::this_mpi_process(communicator))
+           {
+             if (across)
+               {
+                 out << "[Proc" << i << " " << istart << "-" << iend-1 << "]" << ' ';
+                 for (PetscInt i=0; i<nlocal; ++i)
+                   out << val[i] << ' ';
+               }
+             else
+               {
+                 out << "[Proc " << i << " " << istart << "-" << iend-1 << "]" << std::endl;
+                 for (PetscInt i=0; i<nlocal; ++i)
+                   out << val[i] << std::endl;
+               }
+             out << std::endl;
+           }
+       }
+      // reset output format
+      out.flags (old_flags);
+      out.precision(old_precision);
+
+      // restore the representation of the
+      // vector
+      ierr = VecRestoreArray (vector, &val);
+      AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+      AssertThrow (out, ExcIO());
+    }
+
   }
 
 }

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.