]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix Symengine test 8602/head
authorDaniel Arndt <arndtd@ornl.gov>
Mon, 19 Aug 2019 22:53:47 +0000 (18:53 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Wed, 21 Aug 2019 18:11:29 +0000 (14:11 -0400)
include/deal.II/base/tensor.h
tests/symengine/symengine_scalar_operations_04.cc [new file with mode: 0644]
tests/symengine/symengine_scalar_operations_04.output [new file with mode: 0644]

index 23cac8e6459d3badfb8c0c86c80ddaba856d14fc..757c34a473ca743f7910c53ccbc980434efc67a9 100644 (file)
@@ -45,6 +45,13 @@ template <int rank_, int dim, typename Number = double>
 class Tensor;
 template <typename Number>
 class Vector;
+namespace Differentiation
+{
+  namespace SD
+  {
+    class Expression;
+  }
+} // namespace Differentiation
 #endif
 
 #ifndef DOXYGEN
@@ -1341,7 +1348,8 @@ namespace internal
               typename OtherNumber,
               typename std::enable_if<
                 !std::is_integral<
-                  typename ProductType<Number, OtherNumber>::type>::value,
+                  typename ProductType<Number, OtherNumber>::type>::value &&
+                  !std::is_same<Number, Differentiation::SD::Expression>::value,
                 int>::type = 0>
     DEAL_II_CONSTEXPR DEAL_II_CUDA_HOST_DEV inline DEAL_II_ALWAYS_INLINE void
                       division_operator(Tensor<rank, dim, Number> (&t)[dim],
@@ -1360,7 +1368,8 @@ namespace internal
               typename OtherNumber,
               typename std::enable_if<
                 std::is_integral<
-                  typename ProductType<Number, OtherNumber>::type>::value,
+                  typename ProductType<Number, OtherNumber>::type>::value ||
+                  std::is_same<Number, Differentiation::SD::Expression>::value,
                 int>::type = 0>
     DEAL_II_CONSTEXPR DEAL_II_CUDA_HOST_DEV inline DEAL_II_ALWAYS_INLINE void
                       division_operator(dealii::Tensor<rank, dim, Number> (&t)[dim],
diff --git a/tests/symengine/symengine_scalar_operations_04.cc b/tests/symengine/symengine_scalar_operations_04.cc
new file mode 100644 (file)
index 0000000..48826e5
--- /dev/null
@@ -0,0 +1,83 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check that the std::complex overload for the division and multiplication
+// operator works.
+
+#include <deal.II/differentiation/sd.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+namespace SD = Differentiation::SD;
+
+template <typename Number1, typename Number2>
+void
+test()
+{
+  using SD_number_t = SD::Expression;
+
+  const Number1 a = 1.0;
+  const Number2 b = 2.0;
+
+  const SD_number_t symb_a = SD::make_symbol("a");
+  const SD_number_t symb_b = SD::make_symbol("b");
+
+  {
+    const SD_number_t symb_a_division_b = symb_a / b;
+    deallog << "symbolic a/b: " << symb_a_division_b << std::endl;
+
+    const SD_number_t a_division_symb_b = a / symb_b;
+    ;
+    deallog << "a/symbolic b: " << a_division_symb_b << std::endl;
+
+    const SD_number_t substitute_symb_a_division_b =
+      SD::substitute(symb_a_division_b, symb_a, a);
+    deallog << "a/b: " << substitute_symb_a_division_b << std::endl;
+
+    const SD_number_t a_division_substitute_symb_b =
+      SD::substitute(a_division_symb_b, symb_b, b);
+    deallog << "a/b: " << a_division_substitute_symb_b << std::endl;
+  }
+
+  {
+    const SD_number_t symb_a_product_b = symb_a * b;
+    deallog << "symbolic a*b: " << symb_a_product_b << std::endl;
+
+    const SD_number_t a_product_symb_b = a * symb_b;
+    ;
+    deallog << "a*symbolic b: " << a_product_symb_b << std::endl;
+
+    const SD_number_t substitute_symb_a_product_b =
+      SD::substitute(symb_a_product_b, symb_a, a);
+    deallog << "a*b: " << substitute_symb_a_product_b << std::endl;
+
+    const SD_number_t a_product_substitute_symb_b =
+      SD::substitute(a_product_symb_b, symb_b, b);
+    deallog << "a*b: " << a_product_substitute_symb_b << std::endl;
+  }
+}
+
+int
+main()
+{
+  initlog();
+  test<std::complex<float>, std::complex<double>>();
+  test<std::complex<double>, std::complex<float>>();
+  test<std::complex<double>, std::complex<double>>();
+  test<std::complex<float>, std::complex<float>>();
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/symengine/symengine_scalar_operations_04.output b/tests/symengine/symengine_scalar_operations_04.output
new file mode 100644 (file)
index 0000000..d06f7d7
--- /dev/null
@@ -0,0 +1,34 @@
+
+DEAL::symbolic a/b: (0.5 + 0.0*I)*a
+DEAL::a/symbolic b: (1.0 + 0.0*I)/b
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::symbolic a*b: (2.0 + 0.0*I)*a
+DEAL::a*symbolic b: (1.0 + 0.0*I)*b
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::symbolic a/b: (0.5 + 0.0*I)*a
+DEAL::a/symbolic b: (1.0 + 0.0*I)/b
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::symbolic a*b: (2.0 + 0.0*I)*a
+DEAL::a*symbolic b: (1.0 + 0.0*I)*b
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::symbolic a/b: (0.5 + 0.0*I)*a
+DEAL::a/symbolic b: (1.0 + 0.0*I)/b
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::symbolic a*b: (2.0 + 0.0*I)*a
+DEAL::a*symbolic b: (1.0 + 0.0*I)*b
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::symbolic a/b: (0.5 + 0.0*I)*a
+DEAL::a/symbolic b: (1.0 + 0.0*I)/b
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::a/b: 0.5 + 0.0*I
+DEAL::symbolic a*b: (2.0 + 0.0*I)*a
+DEAL::a*symbolic b: (1.0 + 0.0*I)*b
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::a*b: 2.0 + 0.0*I
+DEAL::OK

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.