]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Adress Wolfgang's comments a.k.a add more documentation and switches
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sun, 7 Oct 2018 02:37:41 +0000 (04:37 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 20 Feb 2019 12:46:19 +0000 (13:46 +0100)
include/deal.II/base/utilities.h
include/deal.II/lac/read_write_vector.templates.h
include/deal.II/lac/trilinos_sparse_matrix.h
source/base/utilities.cc

index 71d05faf24ea5beb844768ab26664118bd38921d..dc7414587e540b18a66554067f81aad2b5c3c926 100644 (file)
@@ -888,6 +888,14 @@ namespace Utilities
     const Epetra_Comm &
     comm_self();
 
+    /**
+     * Return a Teuchos::Comm object needed for creation of Tpetra::Maps.
+     *
+     * If deal.II has been configured to use a compiler that does not support
+     * MPI then the resulting communicator will be a serial one. Otherwise,
+     * the communicator will correspond to MPI_COMM_SELF, i.e. a communicator
+     * that comprises only this one processor.
+     */
     const Teuchos::RCP<const Teuchos::Comm<int>> &
     tpetra_comm_self();
 
index 87139aedd776fe82d5a5a013cb7b317dbf90bf2f..9da3c32e68e0218dfa37f42318c17b388a51323c 100644 (file)
@@ -565,30 +565,33 @@ namespace LinearAlgebra
     Assert(size == 0 || values != nullptr, ExcInternalError("Export failed."));
     AssertDimension(size, stored_elements.n_elements());
 
-    if (operation == VectorOperation::insert)
-      {
-        for (size_type i = 0; i < size; ++i)
-          values[i] = new_values[i];
-      }
-    else if (operation == VectorOperation::add)
-      {
-        for (size_type i = 0; i < size; ++i)
-          values[i] += new_values[i];
-      }
-    else if (operation == VectorOperation::min)
+    switch (operation)
       {
-        for (size_type i = 0; i < size; ++i)
-          if (std::real(new_values[i]) - std::real(values[i]) < 0.0)
-            values[i] = new_values[i];
-      }
-    else if (operation == VectorOperation::max)
-      {
-        for (size_type i = 0; i < size; ++i)
-          if (std::real(new_values[i]) - std::real(values[i]) > 0.0)
+        case VectorOperation::insert:
+          for (size_type i = 0; i < size; ++i)
             values[i] = new_values[i];
+          break;
+
+        case VectorOperation::add:
+          for (size_type i = 0; i < size; ++i)
+            values[i] += new_values[i];
+          break;
+
+        case VectorOperation::min:
+          for (size_type i = 0; i < size; ++i)
+            if (std::real(new_values[i]) - std::real(values[i]) < 0.0)
+              values[i] = new_values[i];
+          break;
+
+        case VectorOperation::max:
+          for (size_type i = 0; i < size; ++i)
+            if (std::real(new_values[i]) - std::real(values[i]) > 0.0)
+              values[i] = new_values[i];
+          break;
+
+        default:
+          AssertThrow(false, ExcNotImplemented());
       }
-    else
-      AssertThrow(false, ExcNotImplemented());
   }
 
 
index 2ae42c5ef1576dfa1936ba1252a79d96750aa79f..368fc70a12d4cc560d0cbcd80f520130b6e4cf90 100644 (file)
@@ -1576,6 +1576,9 @@ namespace TrilinosWrappers
      * initialized with the same IndexSet that was used for the column indices
      * of the matrix.
      *
+     * This function will be called when the underlying number type for the
+     * matrix object and the one one for the vector object are the same.
+     *
      * In case of a serial vector, this function will only work when
      * running on one processor, since the matrix object is inherently
      * distributed. Otherwise, an exception will be thrown.
@@ -1585,6 +1588,11 @@ namespace TrilinosWrappers
                                          TrilinosScalar>::value>::type
     vmult(VectorType &dst, const VectorType &src) const;
 
+    /**
+     * Same as the function above for the case that the underlying number type
+     * for the matrix object and the one for the vector object do not coincide.
+     * This case is not implemented. Calling it will result in an runtime error.
+     */
     template <typename VectorType>
     typename std::enable_if<!std::is_same<typename VectorType::value_type,
                                           TrilinosScalar>::value>::type
@@ -1599,12 +1607,20 @@ namespace TrilinosWrappers
      *
      * This function can be called with several types of vector objects,
      * see the discussion about @p VectorType in vmult().
+     *
+     * This function will be called when the underlying number type for the
+     * matrix object and the one one for the vector object are the same.
      */
     template <typename VectorType>
     typename std::enable_if<std::is_same<typename VectorType::value_type,
                                          TrilinosScalar>::value>::type
     Tvmult(VectorType &dst, const VectorType &src) const;
 
+    /**
+     * Same as the function above for the case that the underlying number type
+     * for the matrix object and the one for the vector object do not coincide.
+     * This case is not implemented. Calling it will result in an runtime error.
+     */
     template <typename VectorType>
     typename std::enable_if<!std::is_same<typename VectorType::value_type,
                                           TrilinosScalar>::value>::type
index ff4d0d9573224027df4a67dc5e32826d9846fb79..d5f120192bf6c2002ededdf00eb18a192eba1d5a 100644 (file)
@@ -993,6 +993,8 @@ namespace Utilities
       return communicator;
     }
 
+
+
     const Epetra_Comm &
     comm_self()
     {

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.