]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Tpetra: Fix compiling with complex values 16679/head
authorDaniel Arndt <arndtd@ornl.gov>
Tue, 20 Feb 2024 14:30:15 +0000 (09:30 -0500)
committerDaniel Arndt <arndtd@ornl.gov>
Tue, 20 Feb 2024 14:32:57 +0000 (09:32 -0500)
include/deal.II/lac/read_write_vector.templates.h
include/deal.II/lac/trilinos_tpetra_vector.templates.h

index c72a04c8c7d603d5a6ced2e160fd2493db702ce0..f5ec41ac660d1931b2fc8131c404a054bd0e8695 100644 (file)
@@ -593,12 +593,12 @@ namespace LinearAlgebra
       {
         case VectorOperation::insert:
           for (size_type i = 0; i < size; ++i)
-            values[i] = new_values(i);
+            values[i] = Number(new_values(i));
           break;
 
         case VectorOperation::add:
           for (size_type i = 0; i < size; ++i)
-            values[i] += new_values(i);
+            values[i] += Number(new_values(i));
           break;
 
         case VectorOperation::min:
@@ -609,14 +609,14 @@ namespace LinearAlgebra
           for (size_type i = 0; i < size; ++i)
             {
               Assert(
-                std::imag(new_values(i)) == 0.,
+                std::imag(Number(new_values(i))) == 0.,
                 ExcMessage(
                   "VectorOperation::min is not defined if there is an imaginary part!)"));
               Assert(
                 std::imag(values[i]) == 0.,
                 ExcMessage(
                   "VectorOperation::min is not defined if there is an imaginary part!)"));
-              if (std::real(new_values(i)) - std::real(values[i]) < 0.0)
+              if (std::real(Number(new_values(i))) - std::real(values[i]) < 0.0)
                 values[i] = new_values(i);
             }
           break;
@@ -625,15 +625,15 @@ namespace LinearAlgebra
           for (size_type i = 0; i < size; ++i)
             {
               Assert(
-                std::imag(new_values(i)) == 0.,
+                std::imag(Number(new_values(i))) == 0.,
                 ExcMessage(
                   "VectorOperation::max is not defined if there is an imaginary part!)"));
               Assert(
                 std::imag(values[i]) == 0.,
                 ExcMessage(
                   "VectorOperation::max is not defined if there is an imaginary part!)"));
-              if (std::real(new_values(i)) - std::real(values[i]) > 0.0)
-                values[i] = new_values(i);
+              if (std::real(Number(new_values(i))) - std::real(values[i]) > 0.0)
+                values[i] = Number(new_values(i));
             }
           break;
 
index 04fd6204781f0d138ca4727aef8252eb96f262f1..548c36e1a8d660789830272b108b093c96d7dd0d 100644 (file)
@@ -835,26 +835,34 @@ namespace LinearAlgebra
     bool
     Vector<Number, MemorySpace>::is_non_negative() const
     {
-      // get a representation of the vector and
-      // loop over all the elements
-      Teuchos::ArrayRCP<const Number> data       = vector->getData();
-      const size_type                 n_elements = vector->getLocalLength();
-      unsigned int                    flag       = 0;
-      for (size_type i = 0; i < n_elements; ++i)
+      if constexpr (!std::is_same_v<Number, std::complex<double>> &&
+                    !std::is_same_v<Number, std::complex<float>>)
         {
-          if (data[i] < Number(0))
+          // get a representation of the vector and
+          // loop over all the elements
+          Teuchos::ArrayRCP<const Number> data       = vector->getData();
+          const size_type                 n_elements = vector->getLocalLength();
+          unsigned int                    flag       = 0;
+          for (size_type i = 0; i < n_elements; ++i)
             {
-              flag = 1;
-              break;
+              if (data[i] < Number(0))
+                {
+                  flag = 1;
+                  break;
+                }
             }
-        }
 
-      // Check that the vector is non-negative on _all_ processors.
-      unsigned int num_negative =
-        Utilities::MPI::sum(flag,
-                            Utilities::Trilinos::teuchos_comm_to_mpi_comm(
-                              vector->getMap()->getComm()));
-      return num_negative == 0;
+          // Check that the vector is non-negative on _all_ processors.
+          unsigned int num_negative =
+            Utilities::MPI::sum(flag,
+                                Utilities::Trilinos::teuchos_comm_to_mpi_comm(
+                                  vector->getMap()->getComm()));
+          return num_negative == 0;
+        }
+      Assert(false,
+             ExcMessage("You can't ask a complex value "
+                        "whether it is non-negative."));
+      return true;
     }
 
 

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.