# include <deal.II/differentiation/sd/symengine_number_traits.h>
# include <deal.II/differentiation/sd/symengine_number_types.h>
# include <deal.II/differentiation/sd/symengine_types.h>
+# include <deal.II/differentiation/sd/symengine_utilities.h>
DEAL_II_NAMESPACE_OPEN
convert_expression_vector_to_basic_vector(
const SD::types::symbol_vector &symbol_vector);
+ /**
+ * Extract the symbols (key entries) from a substitution map.
+ *
+ * @note It is guarenteed that the order of extraction of data into the
+ * output vector is the same as that for extract_values().
+ * That is to say that the unzipped key and value pairs as given by
+ * extract_symbols() and extract_values() always have a 1:1
+ * correspondence.
+ */
+ SD::types::symbol_vector
+ extract_symbols(const SD::types::substitution_map &substitution_values);
+
+ /**
+ * Extract the values from a substitution map.
+ * The value entries will be converted into the @p NumberType given
+ * as a template parameter to this function via the @p ExpressionType.
+ *
+ * @note It is guarenteed that the order of extraction of data into the
+ * output vector is the same as that for extract_symbols().
+ * That is to say that the unzipped key and value pairs as given by
+ * extract_symbols() and extract_values() always have a 1:1
+ * correspondence.
+ */
+ template <typename NumberType, typename ExpressionType = SD::Expression>
+ std::vector<NumberType>
+ extract_values(const SD::types::substitution_map &substitution_values);
+
+ /**
+ * Print the key and value pairs stored in a substitution map.
+ */
+ template <typename StreamType>
+ StreamType &
+ print_substitution_map(
+ StreamType & stream,
+ const SD::types::substitution_map &symbol_value_map);
+
} // namespace Utilities
} // namespace SD
} // namespace Differentiation
+/* -------------------- inline and template functions ------------------ */
+
+
+# ifndef DOXYGEN
+
+
+namespace Differentiation
+{
+ namespace SD
+ {
+ namespace Utilities
+ {
+ template <typename NumberType, typename ExpressionType>
+ std::vector<NumberType>
+ extract_values(const SD::types::substitution_map &substitution_values)
+ {
+ std::vector<NumberType> values;
+ values.reserve(substitution_values.size());
+
+ for (const auto &entry : substitution_values)
+ values.push_back(
+ static_cast<NumberType>(ExpressionType(entry.second)));
+
+ return values;
+ }
+
+
+ template <typename StreamType>
+ StreamType &
+ print_substitution_map(
+ StreamType & stream,
+ const SD::types::substitution_map &symbol_value_map)
+ {
+ for (const auto &entry : symbol_value_map)
+ stream << entry.first << " = " << entry.second << "\n";
+
+ stream << std::flush;
+ return stream;
+ }
+
+ } // namespace Utilities
+
+ } // namespace SD
+} // namespace Differentiation
+
+
+# endif // DOXYGEN
+
+
DEAL_II_NAMESPACE_CLOSE
#endif // DEAL_II_WITH_SYMENGINE
return symb_vec;
}
+
+ SD::types::symbol_vector
+ extract_symbols(const SD::types::substitution_map &substitution_values)
+ {
+ SD::types::symbol_vector symbols;
+ symbols.reserve(substitution_values.size());
+
+ for (typename SD::types::substitution_map::const_iterator it =
+ substitution_values.begin();
+ it != substitution_values.end();
+ ++it)
+ symbols.push_back(it->first);
+
+ return symbols;
+ }
+
} // namespace Utilities
} // namespace SD
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 substitution map decomposition and printing work as expected
+
+#include <deal.II/differentiation/sd.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+namespace SD = Differentiation::SD;
+
+
+int
+main()
+{
+ initlog();
+
+ using SD_number_t = SD::Expression;
+
+ const SD::types::substitution_map sub_map{
+ {SD_number_t("c"), SD_number_t(1.0)},
+ {SD_number_t("b"), SD_number_t(2)},
+ {SD_number_t("a"), SD_number_t(3.0f)}};
+
+ const SD::types::symbol_vector symbols =
+ SD::Utilities::extract_symbols(sub_map);
+ const std::vector<double> values =
+ SD::Utilities::extract_values<double>(sub_map);
+ Assert(values.size() == symbols.size(),
+ ExcDimensionMismatch(values.size(), symbols.size()));
+
+ // Print the map itself (this should be dictionary ordered)
+ deallog << "Print substitution map" << std::endl;
+ SD::Utilities::print_substitution_map(deallog, sub_map);
+
+ // Print the extracted symbol-value pairs (the ordering should match that of
+ // the map)
+ deallog << "Print extracted symbol-value pairs" << std::endl;
+ for (unsigned int i = 0; i < symbols.size(); ++i)
+ deallog << symbols[i] << " = " << values[i] << std::endl;
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::Print substitution map
+DEAL::a = 3.0
+b = 2
+c = 1.0
+Print extracted symbol-value pairs
+DEAL::a = 3.00000
+DEAL::b = 2.00000
+DEAL::c = 1.00000
+DEAL::OK