]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix issues discovered by Coverity 5854/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 3 Feb 2018 21:08:13 +0000 (22:08 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sun, 4 Feb 2018 19:05:58 +0000 (20:05 +0100)
examples/step-32/step-32.cc
examples/step-50/step-50.cc
include/deal.II/lac/la_parallel_vector.templates.h
source/distributed/tria.cc
source/dofs/dof_handler_policy.cc
source/sundials/arkode.cc
source/sundials/ida.cc
source/sundials/kinsol.cc

index a8d6c7019b264de90aedf7fd8b9a734d55aafcb4..bc88b31738ada12a26569f08f35210efd906fdb4 100644 (file)
@@ -3726,14 +3726,14 @@ start_time_iteration:
 // by changing the constant dimension below to 3.
 int main (int argc, char *argv[])
 {
-  using namespace Step32;
-  using namespace dealii;
-
-  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv,
-                                                      numbers::invalid_unsigned_int);
-
   try
     {
+      using namespace Step32;
+      using namespace dealii;
+
+      Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv,
+                                                          numbers::invalid_unsigned_int);
+
       std::string parameter_filename;
       if (argc>=2)
         parameter_filename = argv[1];
index 6bff73413cc19813318379a9fbd6fb3d0170889a..22258225fc6478729f81a20265dd369da8e66fd7 100644 (file)
@@ -968,13 +968,13 @@ namespace Step50
 // in step-6:
 int main (int argc, char *argv[])
 {
-  dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
-
   try
     {
       using namespace dealii;
       using namespace Step50;
 
+      Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
       LaplaceProblem<2> laplace_problem(1/*degree*/);
       laplace_problem.run ();
     }
index 0ab50f602cfc25c28f90ffc4814816c6face6fbf..6cb5794d83a2f888b5e71627e55df4752033a02b 100644 (file)
@@ -290,7 +290,12 @@ namespace LinearAlgebra
     inline
     Vector<Number>::~Vector ()
     {
-      clear_mpi_requests();
+      try
+        {
+          clear_mpi_requests();
+        }
+      catch (...)
+        {}
     }
 
 
index 95bfc1a108b65e85582cc105046d3ebdf90fbd1c..4195b3673b4bba1f199cf53606f4d794da16c9eb 100644 (file)
@@ -1307,12 +1307,17 @@ namespace parallel
       // virtual functions called in constructors and destructors never use the
       // override in a derived class
       // for clarity be explicit on which function is called
-      Triangulation<dim, spacedim>::clear ();
+      try
+        {
+          Triangulation<dim, spacedim>::clear ();
+        }
+      catch (...)
+        {}
 
-      Assert (triangulation_has_content == false,
-              ExcInternalError());
-      Assert (connectivity == nullptr,    ExcInternalError());
-      Assert (parallel_forest == nullptr, ExcInternalError());
+      AssertNothrow (triangulation_has_content == false,
+                     ExcInternalError());
+      AssertNothrow (connectivity == nullptr,    ExcInternalError());
+      AssertNothrow (parallel_forest == nullptr, ExcInternalError());
     }
 
 
@@ -3741,7 +3746,7 @@ namespace parallel
     template <int spacedim>
     Triangulation<1,spacedim>::~Triangulation ()
     {
-      Assert (false, ExcNotImplemented());
+      AssertNothrow (false, ExcNotImplemented());
     }
 
 
index 03c9bae9e3e730969a161e489e5b9dfeb09908a0..a56e9a12b45b1dbb53c9e5a4955d718321587cb3 100644 (file)
@@ -3638,11 +3638,17 @@ namespace internal
           // different tags for phase 1 and 2, but the cost of a
           // barrier is negligible compared to everything else we do
           // here
-          const parallel::distributed::Triangulation< dim, spacedim > *triangulation
-            = (dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim>*>
-               (&dof_handler.get_triangulation()));
-          const int ierr = MPI_Barrier(triangulation->get_communicator());
-          AssertThrowMPI(ierr);
+          if (const auto *triangulation =
+                dynamic_cast<const parallel::distributed::Triangulation<dim,spacedim>*>(&dof_handler.get_triangulation()))
+            {
+              const int ierr = MPI_Barrier(triangulation->get_communicator());
+              AssertThrowMPI(ierr);
+            }
+          else
+            {
+              Assert(false, ExcMessage("The function communicate_dof_indices_on_marked_cells() "
+                                       "only works with parallel distributed triangulations."));
+            }
 #endif
         }
 
index 68ebed5cb8405b3eb9dad4e7dbba2df71cb03ca8..205c0d9b474656dd86adb018dd100753a257d24f 100644 (file)
@@ -221,6 +221,8 @@ namespace SUNDIALS
                              const MPI_Comm mpi_comm) :
     data(data),
     arkode_mem(nullptr),
+    yy(nullptr),
+    abs_tolls(nullptr),
     communicator(is_serial_vector<VectorType>::value ?
                  MPI_COMM_SELF :
                  Utilities::MPI::duplicate_communicator(mpi_comm))
@@ -237,7 +239,7 @@ namespace SUNDIALS
     if (is_serial_vector<VectorType>::value == false)
       {
         const int ierr = MPI_Comm_free(&communicator);
-        AssertThrowMPI(ierr);
+        AssertNothrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
       }
 #endif
   }
index 160bb1e3a1dc31bdddf8c8fdc0e9d091e624ce9d..11afc8fc8841e5a255764a411c86c67204c28ae8 100644 (file)
@@ -147,6 +147,10 @@ namespace SUNDIALS
                        const MPI_Comm mpi_comm) :
     data(data),
     ida_mem(nullptr),
+    yy(nullptr),
+    yp(nullptr),
+    abs_tolls(nullptr),
+    diff_id(nullptr),
     communicator(is_serial_vector<VectorType>::value ?
                  MPI_COMM_SELF :
                  Utilities::MPI::duplicate_communicator(mpi_comm))
@@ -163,7 +167,7 @@ namespace SUNDIALS
     if (is_serial_vector<VectorType>::value == false)
       {
         const int ierr = MPI_Comm_free(&communicator);
-        AssertThrowMPI(ierr);
+        AssertNothrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
       }
 #endif
   }
index fe00b17db2cbf3e7a4815ed1a3bc918fd8ec8488..ba07fb638f69e35f7048da8503147fbea73f655e 100644 (file)
@@ -148,6 +148,9 @@ namespace SUNDIALS
   KINSOL<VectorType>::KINSOL(const AdditionalData &data, const MPI_Comm mpi_comm) :
     data(data),
     kinsol_mem(nullptr),
+    solution(nullptr),
+    u_scale(nullptr),
+    f_scale(nullptr),
     communicator(is_serial_vector<VectorType>::value ?
                  MPI_COMM_SELF :
                  Utilities::MPI::duplicate_communicator(mpi_comm))
@@ -166,7 +169,7 @@ namespace SUNDIALS
     if (is_serial_vector<VectorType>::value == false)
       {
         const int ierr = MPI_Comm_free(&communicator);
-        AssertThrowMPI(ierr);
+        AssertNothrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
       }
 #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.