]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Allow SymmetricGradient::operator* for aggregates
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 28 Jan 2017 11:22:59 +0000 (12:22 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 28 Jan 2017 16:42:50 +0000 (17:42 +0100)
include/deal.II/base/symmetric_tensor.h
tests/base/symmetric_tensor_38.cc [new file with mode: 0644]
tests/base/symmetric_tensor_38.output [new file with mode: 0644]

index 2a5b4c16fae5343741ad2032003876b52b558bfc..65b8626f0a1d482ceab5af2be8a563fbe69b206c 100644 (file)
@@ -2916,7 +2916,9 @@ operator * (const SymmetricTensor<rank,dim,Number> &t,
   // (as well as with switched arguments and double<->float).
   typedef typename ProductType<Number,OtherNumber>::type product_type;
   SymmetricTensor<rank,dim,product_type> tt(t);
-  tt *= product_type(factor);
+  product_type new_factor;
+  new_factor = factor;
+  tt *= new_factor;
   return tt;
 }
 
diff --git a/tests/base/symmetric_tensor_38.cc b/tests/base/symmetric_tensor_38.cc
new file mode 100644 (file)
index 0000000..ad780e8
--- /dev/null
@@ -0,0 +1,75 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 operator* with a VectorizedArray
+
+#include "../tests.h"
+#include <deal.II/base/symmetric_tensor.h>
+#include <deal.II/base/vectorization.h>
+#include <deal.II/base/tensor.h>
+#include <deal.II/base/logstream.h>
+
+#include <fstream>
+#include <iomanip>
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog << std::setprecision(5);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  const int dim = 3;
+  SymmetricTensor<2,dim,VectorizedArray<float> > s1;
+  SymmetricTensor<2,dim,VectorizedArray<double> > s2;
+
+  for (unsigned int i=0; i<dim; ++i)
+    for (unsigned int j=0; j<dim; ++j)
+      {
+        s1[i][j] = 1+j+i*dim;
+        s2[i][j] = 2+(j+i)*dim;
+      }
+
+  float factor_float = 2.f;
+  double factor_double = 3.;
+  SymmetricTensor<2,dim,VectorizedArray<float> > r1 = factor_float*s1;
+  SymmetricTensor<2,dim,VectorizedArray<double> > r2 = factor_double*s2;
+
+  for (unsigned int i=0; i<dim; ++i)
+    for (unsigned int j=0; j<dim; ++j)
+      {
+        deallog << i << "\t" << j << std::endl;
+        deallog << r1[i][j][0] << std::endl;
+        AssertThrow(std::abs(r1[i][j][0]-factor_float*s1[i][j][0])<1.e-10,
+                    ExcInternalError());
+        for (unsigned int k=1; k<VectorizedArray<float>::n_array_elements; ++k)
+          {
+            AssertThrow(std::abs(r1[i][j][k]-r1[i][j][0])<1.e-10,
+                        ExcInternalError());
+          }
+        deallog << r2[i][j][0] << std::endl;
+        Assert(std::abs(r2[i][j][0]-factor_double*s2[i][j][0])<1.e-10,
+               ExcInternalError());
+        for (unsigned int k=1; k<VectorizedArray<double>::n_array_elements; ++k)
+          {
+            AssertThrow(std::abs(r2[i][j][k]-r2[i][j][0])<1.e-10,
+                        ExcInternalError());
+          }
+        deallog << std::endl;
+      }
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/base/symmetric_tensor_38.output b/tests/base/symmetric_tensor_38.output
new file mode 100644 (file)
index 0000000..76d7764
--- /dev/null
@@ -0,0 +1,38 @@
+
+DEAL::0        0
+DEAL::2.0000
+DEAL::6.0000
+DEAL::
+DEAL::0        1
+DEAL::8.0000
+DEAL::15.000
+DEAL::
+DEAL::0        2
+DEAL::14.000
+DEAL::24.000
+DEAL::
+DEAL::1        0
+DEAL::8.0000
+DEAL::15.000
+DEAL::
+DEAL::1        1
+DEAL::10.000
+DEAL::24.000
+DEAL::
+DEAL::1        2
+DEAL::16.000
+DEAL::33.000
+DEAL::
+DEAL::2        0
+DEAL::14.000
+DEAL::24.000
+DEAL::
+DEAL::2        1
+DEAL::16.000
+DEAL::33.000
+DEAL::
+DEAL::2        2
+DEAL::18.000
+DEAL::42.000
+DEAL::
+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.