]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add functions to create SD subsitution maps from tensor types 7973/head
authorJean-Paul Pelteret <jppelteret@gmail.com>
Sat, 27 Apr 2019 09:42:14 +0000 (11:42 +0200)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Sun, 28 Apr 2019 15:22:47 +0000 (17:22 +0200)
include/deal.II/differentiation/sd/symengine_tensor_operations.h
tests/symengine/substitution_maps_tensor_02.cc [new file with mode: 0644]
tests/symengine/substitution_maps_tensor_02.output [new file with mode: 0644]

index c325bf9cc89c412c8b3a0bbe709e939b872fa597..3c148486975b60adcaa0c0009840af68196ea5a5 100644 (file)
@@ -366,6 +366,44 @@ namespace Differentiation
 
     //@}
 
+    /**
+     * @name Symbol substitution map creation
+     */
+    //@{
+
+    /**
+     * Return a substitution map that has the entry keys given by the
+     * @p symbol_tensor and the values given by the @p value_tensor. It is
+     * expected that all key entries be valid symbols or symbolic expressions.
+     *
+     * It is possible to map symbolic types to other symbolic types
+     * using this function. For more details on this, see the other
+     * \ref make_substitution_map(const Expression &,const ValueType &)
+     * function.
+     */
+    template <int rank, int dim, typename ValueType>
+    types::substitution_map
+    make_substitution_map(const Tensor<rank, dim, Expression> &symbol_tensor,
+                          const Tensor<rank, dim, ValueType> & value_tensor);
+
+    /**
+     * Return a substitution map that has the entry keys given by the
+     * @p symbol_tensor and the values given by the @p value_tensor. It is
+     * expected that all key entries be valid symbols or symbolic expressions.
+     *
+     * It is possible to map symbolic types to other symbolic types
+     * using this function. For more details on this, see the other
+     * \ref make_substitution_map(const Expression &,const ValueType &)
+     * function.
+     */
+    template <int rank, int dim, typename ValueType>
+    types::substitution_map
+    make_substitution_map(
+      const SymmetricTensor<rank, dim, Expression> &symbol_tensor,
+      const SymmetricTensor<rank, dim, ValueType> & value_tensor);
+
+    //@}
+
     /**
      * @name Symbol substitution map enlargement
      */
@@ -811,6 +849,32 @@ namespace Differentiation
     }
 
 
+    /* ---------------- Symbolic substitution map creation --------------*/
+
+
+    template <int rank, int dim, typename ValueType>
+    types::substitution_map
+    make_substitution_map(const Tensor<rank, dim, Expression> &symbol_tensor,
+                          const Tensor<rank, dim, ValueType> & value_tensor)
+    {
+      types::substitution_map substitution_map;
+      add_to_substitution_map(substitution_map, symbol_tensor, value_tensor);
+      return substitution_map;
+    }
+
+
+    template <int rank, int dim, typename ValueType>
+    types::substitution_map
+    make_substitution_map(
+      const SymmetricTensor<rank, dim, Expression> &symbol_tensor,
+      const SymmetricTensor<rank, dim, ValueType> & value_tensor)
+    {
+      types::substitution_map substitution_map;
+      add_to_substitution_map(substitution_map, symbol_tensor, value_tensor);
+      return substitution_map;
+    }
+
+
     /* ---------------- Symbolic substitution map enlargement --------------*/
 
 
diff --git a/tests/symengine/substitution_maps_tensor_02.cc b/tests/symengine/substitution_maps_tensor_02.cc
new file mode 100644 (file)
index 0000000..cc44eb8
--- /dev/null
@@ -0,0 +1,61 @@
+// ---------------------------------------------------------------------
+//
+// 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 functions to create substitution maps from tensor variables work
+// correctly.
+
+#include <deal.II/differentiation/sd.h>
+
+#include <string>
+
+#include "../tests.h"
+
+using namespace dealii;
+namespace SD = Differentiation::SD;
+namespace SE = ::SymEngine;
+
+int
+main()
+{
+  initlog();
+  const unsigned int dim = 2;
+
+  Tensor<2, dim>          t;
+  SymmetricTensor<2, dim> st;
+  for (unsigned int i = 0; i < t.n_independent_components; ++i)
+    t[t.unrolled_to_component_indices(i)] = i;
+  for (unsigned int i = 0; i < st.n_independent_components; ++i)
+    st[st.unrolled_to_component_indices(i)] = i;
+
+
+  SD::types::substitution_map substitution_map;
+  SD::add_to_substitution_map(substitution_map,
+                              SD::Expression("x1"),
+                              SD::Expression(1));
+
+  SD::merge_substitution_maps(
+    substitution_map,
+    SD::make_substitution_map(SD::make_tensor_of_symbols<2, dim>("t1"), t),
+    SD::make_substitution_map(
+      SD::make_symmetric_tensor_of_symbols<2, dim>("st1"), st),
+    SD::make_substitution_map(SD::make_tensor_of_symbols<2, dim>("t2"), t),
+    SD::make_substitution_map(
+      SD::make_symmetric_tensor_of_symbols<2, dim>("st2"), st));
+
+  SD::Utilities::print_substitution_map(deallog, substitution_map);
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/symengine/substitution_maps_tensor_02.output b/tests/symengine/substitution_maps_tensor_02.output
new file mode 100644 (file)
index 0000000..c562f76
--- /dev/null
@@ -0,0 +1,17 @@
+
+DEAL::x1 = 1
+t2_01 = 1.0
+t2_00 = 0.0
+t2_10 = 2.0
+t2_11 = 3.0
+t1_00 = 0.0
+t1_01 = 1.0
+t1_11 = 3.0
+t1_10 = 2.0
+st2_11 = 1.0
+st2_00 = 0.0
+st2_01 = 2.0
+st1_11 = 1.0
+st1_00 = 0.0
+st1_01 = 2.0
+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.