]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
MPI_InitFinalize can now also initialize PETSc.
authorheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 3 Nov 2012 03:27:29 +0000 (03:27 +0000)
committerheister <heister@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 3 Nov 2012 03:27:29 +0000 (03:27 +0000)
git-svn-id: https://svn.dealii.org/trunk@27328 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/examples/step-17/step-17.cc
deal.II/examples/step-18/step-18.cc
deal.II/examples/step-40/step-40.cc
deal.II/include/deal.II/base/mpi.h
deal.II/source/base/mpi.cc

index 9ff00210fb48471fea2341ce4af69b9fcbe7efbe..ad8e51d12e50470b80534373028e8d208b85ae02 100644 (file)
@@ -108,6 +108,11 @@ never working correctly and it is not used.
 
 
 <ol>
+<li> The class Utilities::MPI::MPI_InitFinalize now also initializes
+PETSc, when PETSc is installed.
+<br>
+(Timo Heister, 2012/11/02)
+
 <li> step-6 now uses ConstraintMatrix::distribute_local_to_global()
 instead of condense(), which is the preferred way to use a ConstraintMatrix
  (and the only sensible way in parallel).
index a93c3849a830f5cf7bca59a9fb2606c7882b0828..bdb6faaed6c37f08bce3812d7e4fec788b66dd07 100644 (file)
@@ -1226,14 +1226,14 @@ int main (int argc, char **argv)
                                        // Here is the only real difference:
                                        // PETSc requires that we initialize it
                                        // at the beginning of the program, and
-                                       // un-initialize it at the end. So we
-                                       // call <code>PetscInitialize</code> and
-                                       // <code>PetscFinalize</code>. The original code
+                                       // un-initialize it at the end. The
+                                       // class MPI_InitFinalize takes care
+                                       // of that. The original code
                                        // sits in between, enclosed in braces
                                        // to make sure that the
                                        // <code>elastic_problem</code> variable goes
                                        // out of scope (and is destroyed)
-                                       // before we call
+                                       // before PETSc is closed with
                                        // <code>PetscFinalize</code>. (If we wouldn't
                                        // use braces, the destructor of
                                        // <code>elastic_problem</code> would run after
@@ -1241,7 +1241,7 @@ int main (int argc, char **argv)
                                        // destructor involves calls to PETSc
                                        // functions, we would get strange
                                        // error messages from PETSc.)
-      PetscInitialize(&argc,&argv,0,0);
+      Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
 
       {
         deallog.depth_console (0);
@@ -1249,8 +1249,6 @@ int main (int argc, char **argv)
         ElasticProblem<2> elastic_problem;
         elastic_problem.run ();
       }
-
-      PetscFinalize();
     }
   catch (std::exception &exc)
     {
index bc051f79c6dea0dbd1c69d1a436d29cb08d47ead..e7e1cd304b520836167d151fecd4b51afdd40bbe 100644 (file)
@@ -2821,16 +2821,13 @@ int main (int argc, char **argv)
       using namespace dealii;
       using namespace Step18;
 
-      PetscInitialize(&argc,&argv,0,0);
-
+      Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
       {
         deallog.depth_console (0);
 
         TopLevel<3> elastic_problem;
         elastic_problem.run ();
       }
-
-      PetscFinalize();
     }
   catch (std::exception &exc)
     {
index 3b201990e59ac87f63d5a01e3c3e6f457fd31a27..229b562d298b997c623fc4abc847d1447add9d9b 100644 (file)
@@ -881,17 +881,16 @@ namespace Step40
                                  // step-6. Like in the other programs
                                  // that use PETSc, we have to
                                  // inialize and finalize PETSc, which
-                                 // also initializes and finalizes the
-                                 // MPI subsystem.
+                                 // is done using the helper object
+                                 // MPI_InitFinalize.
                                  //
                                  // Note how we enclose the use the
                                  // use of the LaplaceProblem class in
                                  // a pair of braces. This makes sure
                                  // that all member variables of the
                                  // object are destroyed by the time
-                                 // we hit the
-                                 // <code>PetscFinalize</code>
-                                 // call. Not doing this will lead to
+                                 // we destroy the mpi_intialization
+                                 // object. Not doing this will lead to
                                  // strange and hard to debug errors
                                  // when <code>PetscFinalize</code>
                                  // first deletes all PETSc vectors
@@ -906,15 +905,13 @@ int main(int argc, char *argv[])
       using namespace dealii;
       using namespace Step40;
 
-      PetscInitialize(&argc, &argv, PETSC_NULL, PETSC_NULL);
+      Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
       deallog.depth_console (0);
 
       {
         LaplaceProblem<2> laplace_problem_2d;
         laplace_problem_2d.run ();
       }
-
-      PetscFinalize();
     }
   catch (std::exception &exc)
     {
index 2fdb14da05da9a65957558e45f9588e0b05fb6fa..6cf6a981ae2dca6e3b9d02b1d6843eeb2ce00886 100644 (file)
@@ -267,6 +267,12 @@ namespace Utilities
                                       * program and to shut it down again at
                                       * the end.
                                       *
+                                      * If deal.II is configured with PETSc,
+                                      * the library will also be initialized
+                                      * in the beginning and destructed at the
+                                      * end automatically (internally by calling
+                                      * PetscInitialize() and PetscFinalize()).
+                                      *
                                       * If a program uses MPI one would
                                       * typically just create an object of
                                       * this type at the beginning of
index fe721b947e3b159a15f812366c0a831c467902d4..c72910dc308945ed0686ff42fed82fcd4ecede6a 100644 (file)
 #  endif
 #endif
 
+#ifdef DEAL_II_USE_PETSC
+#  ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+#    include <petscsys.h>
+#  endif
+#endif
+
+
 DEAL_II_NAMESPACE_OPEN
 
 
@@ -299,6 +306,11 @@ namespace Utilities
                           "in a program since it initializes the MPI system."));
 
 #ifdef DEAL_II_COMPILER_SUPPORTS_MPI
+      // if we have PETSc, we will initialize it and let it handle MPI.
+      // Otherwise, we will do it.
+#ifdef DEAL_II_USE_PETSC
+      PetscInitialize(&argc, &argv, PETSC_NULL, PETSC_NULL);
+#else
       int MPI_has_been_started = 0;
       MPI_Initialized(&MPI_has_been_started);
       AssertThrow (MPI_has_been_started == 0,
@@ -308,6 +320,7 @@ namespace Utilities
       mpi_err = MPI_Init (&argc, &argv);
       AssertThrow (mpi_err == 0,
                    ExcMessage ("MPI could not be initialized."));
+#endif
 #else
                                        // make sure the compiler doesn't warn
                                        // about these variables
@@ -350,6 +363,11 @@ namespace Utilities
         ::release_unused_memory ();
 #  endif
 
+#ifdef DEAL_II_USE_PETSC
+      PetscFinalize();
+#else
+
+
       int mpi_err = 0;
 
       int MPI_has_been_started = 0;
@@ -371,6 +389,7 @@ namespace Utilities
 
       AssertThrow (mpi_err == 0,
                    ExcMessage ("An error occurred while calling MPI_Finalize()"));
+#endif
 #endif
     }
 

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.