From: Roland Date: Fri, 7 Jun 2019 11:35:35 +0000 (+0200) Subject: Elementwise multiplication for Tensor X-Git-Tag: v9.2.0-rc1~1392^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=3a1eb80b2d5550df2eeb341a37c6f216a82e648a;p=dealii.git Elementwise multiplication for Tensor --- diff --git a/doc/news/changes/minor/20190610RolandRichter b/doc/news/changes/minor/20190610RolandRichter new file mode 100644 index 0000000000..47913c89ff --- /dev/null +++ b/doc/news/changes/minor/20190610RolandRichter @@ -0,0 +1,4 @@ +Added: The function "Schur-product" was added for the Tensor-class, +allowing the entrywise multiplication of tensor objects of general rank +
+(Roland Richter, 2019/06/10) diff --git a/include/deal.II/base/tensor.h b/include/deal.II/base/tensor.h index 140b2e1df4..2531af9fad 100644 --- a/include/deal.II/base/tensor.h +++ b/include/deal.II/base/tensor.h @@ -1817,6 +1817,73 @@ DEAL_II_CONSTEXPR DEAL_II_CUDA_HOST_DEV inline DEAL_II_ALWAYS_INLINE return tmp; } +/** + * Entrywise multiplication of two tensor objects of rank 0 (i.e. the + * multiplication of two scalar values). + * + * @relatesalso Tensor + */ +template +inline DEAL_II_ALWAYS_INLINE + Tensor<0, dim, typename ProductType::type> + schur_product(const Tensor<0, dim, Number> & src1, + const Tensor<0, dim, OtherNumber> &src2) +{ + Tensor<0, dim, typename ProductType::type> tmp(src1); + + tmp *= src2; + + return tmp; +} + +/** + * Entrywise multiplication of two tensor objects of rank 1. + * + * @relatesalso Tensor + */ +template +inline DEAL_II_ALWAYS_INLINE + Tensor<1, dim, typename ProductType::type> + schur_product(const Tensor<1, dim, Number> & src1, + const Tensor<1, dim, OtherNumber> &src2) +{ + Tensor<1, dim, typename ProductType::type> tmp(src1); + + for (unsigned int i = 0; i < dim; ++i) + tmp[i] *= src2[i]; + + return tmp; +} + +/** + * Entrywise multiplication of two tensor objects of general rank. + * + * This multiplication is also called "Hadamard-product" (c.f. + * https://en.wikipedia.org/wiki/Hadamard_product_(matrices)), and generates a + * new tensor of size : + * @f[ + * \text{result}_{i, j} + * = \text{left}_{i, j}\cdot + * \text{right}_{i, j} + * @f] + * + * @tparam rank The rank of both tensors. + * + * @relatesalso Tensor + */ +template +inline DEAL_II_ALWAYS_INLINE + Tensor::type> + schur_product(const Tensor & src1, + const Tensor &src2) +{ + Tensor::type> tmp(src1); + + for (unsigned int i = 0; i < dim; ++i) + tmp[i] = schur_product(src1[i], src2[i]); + + return tmp; +} //@} /** diff --git a/tests/tensors/tensor_10.cc b/tests/tensors/tensor_10.cc new file mode 100644 index 0000000000..93b5085a76 --- /dev/null +++ b/tests/tensors/tensor_10.cc @@ -0,0 +1,98 @@ +// --------------------------------------------------------------------- +// +// 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.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + +// check schur_product for Tensor with rank = 0, 1, 2, 3 and dim = 2 + +#include + +#include + +#include "../tests.h" + +template +void +test_same(const Tensor &t, const Tensor &compare) +{ + deallog << "Constant rank " << rank << " and dim " << dim << '\t' + << schur_product(t, t) << " compare " << compare << std::endl; +} + +template +void +test_different(const Tensor &t, + const Tensor &t2, + const Tensor &compare) +{ + deallog << "Constant rank " << rank << " and dim " << dim << '\t' + << schur_product(t, t2) << " compare " << compare << std::endl; +} + +int +main() +{ + initlog(); + + { + Tensor<0, 2> a, b, compare; + a = 2; + b = 2; + compare = 4; + test_same(a, compare); + test_different(a, b, compare); + } + + { + Tensor<1, 2> a, b, compare; + for (size_t d = 0; d < 2; ++d) + { + a[d] = 2; + b[d] = 2; + compare[d] = 4; + } + test_same(a, compare); + test_different(a, b, compare); + } + + { + Tensor<2, 2> a, b, compare; + for (size_t c = 0; c < 2; ++c) + for (size_t d = 0; d < 2; ++d) + { + a[c][d] = 2; + b[c][d] = 2; + compare[c][d] = 4; + } + test_same(a, compare); + test_different(a, b, compare); + } + + { + Tensor<3, 2> a, b, compare; + + for (size_t c = 0; c < 2; ++c) + for (size_t d = 0; d < 2; ++d) + for (size_t e = 0; e < 2; ++e) + { + a[c][d][e] = 2; + b[c][d][e] = 2; + compare[c][d][e] = 4; + } + test_same(a, compare); + test_different(a, b, compare); + } + + return 0; +} diff --git a/tests/tensors/tensor_10.output b/tests/tensors/tensor_10.output new file mode 100644 index 0000000000..baa36bcc51 --- /dev/null +++ b/tests/tensors/tensor_10.output @@ -0,0 +1,9 @@ + +DEAL::Constant rank 0 and dim 2 4.00000 compare 4.00000 +DEAL::Constant rank 0 and dim 2 4.00000 compare 4.00000 +DEAL::Constant rank 1 and dim 2 4.00000 4.00000 compare 4.00000 4.00000 +DEAL::Constant rank 1 and dim 2 4.00000 4.00000 compare 4.00000 4.00000 +DEAL::Constant rank 2 and dim 2 4.00000 4.00000 4.00000 4.00000 compare 4.00000 4.00000 4.00000 4.00000 +DEAL::Constant rank 2 and dim 2 4.00000 4.00000 4.00000 4.00000 compare 4.00000 4.00000 4.00000 4.00000 +DEAL::Constant rank 3 and dim 2 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 compare 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 +DEAL::Constant rank 3 and dim 2 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 compare 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000 4.00000