--- /dev/null
+## ------------------------------------------------------------------------
+##
+## SPDX-License-Identifier: LGPL-2.1-or-later
+## Copyright (C) 2012 - 2022 by the deal.II authors
+##
+## This file is part of the deal.II library.
+##
+## Part of the source code is dual licensed under Apache-2.0 WITH
+## LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+## governing the source code and code contributions can be found in
+## LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+##
+## ------------------------------------------------------------------------
+
+#
+# Configuration for the MAGIC_ENUM library:
+#
+
+configure_feature(MAGIC_ENUM)
--- /dev/null
+## ------------------------------------------------------------------------
+##
+## SPDX-License-Identifier: LGPL-2.1-or-later
+## Copyright (C) 2012 - 2022 by the deal.II authors
+##
+## This file is part of the deal.II library.
+##
+## Part of the source code is dual licensed under Apache-2.0 WITH
+## LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+## governing the source code and code contributions can be found in
+## LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+##
+## ------------------------------------------------------------------------
+
+#
+# Try to find the MAGIC ENUM library
+#
+# This module exports
+#
+# MAGIC_ENUM_FOUND
+# MAGIC_ENUM_INCLUDE_DIRS
+#
+
+set(MAGIC_ENUM_DIR "" CACHE PATH "An optional hint to a MAGIC_ENUM installation")
+set_if_empty(MAGIC_ENUM_DIR "$ENV{MAGIC_ENUM_DIR}")
+
+deal_ii_find_path(MAGIC_ENUM_INCLUDE_DIR magic_enum.hpp
+ HINTS ${MAGIC_ENUM_DIR}
+ PATH_SUFFIXES include
+ )
+
+process_feature(MAGIC_ENUM
+ LIBRARIES
+ INCLUDE_DIRS
+ REQUIRED MAGIC_ENUM_INCLUDE_DIR
+ CLEAR MAGIC_ENUM_INCLUDE_DIR
+ )
--- /dev/null
+New: Added support for magic_enum.hpp. Now PatternsTools::Convert<T> works also when T is an enum type.
+<br>
+(Luca Heltai, 2024/08/12)
#cmakedefine DEAL_II_WITH_LAPACK
#cmakedefine LAPACK_WITH_64BIT_BLAS_INDICES
#cmakedefine DEAL_II_LAPACK_WITH_MKL
+#cmakedefine DEAL_II_WITH_MAGIC_ENUM
#cmakedefine DEAL_II_WITH_METIS
#cmakedefine DEAL_II_WITH_MPI
#cmakedefine DEAL_II_WITH_MUPARSER
#include <utility>
#include <vector>
+#ifdef DEAL_II_WITH_MAGIC_ENUM
+# include <magic_enum.hpp>
+#endif
+
DEAL_II_NAMESPACE_OPEN
// forward declarations for interfaces and friendship
}
};
+#ifdef DEAL_II_WITH_MAGIC_ENUM
+ // Enums
+ template <class T>
+ struct Convert<T, typename std::enable_if<std::is_enum<T>::value>::type>
+ {
+ static std::unique_ptr<Patterns::PatternBase>
+ to_pattern()
+ {
+ const auto n = magic_enum::enum_names<T>();
+ std::vector<std::string> names = {n.begin(), n.end()};
+ const auto selection =
+ Patterns::Tools::Convert<decltype(names)>::to_string(
+ names,
+ Patterns::List(
+ Patterns::Anything(), names.size(), names.size(), "|"));
+ // Allow parsing a list of enums, and make bitwise or between them
+ return Patterns::List(Patterns::Selection(selection),
+ 0,
+ names.size(),
+ "|")
+ .clone();
+ }
+
+ static std::string
+ to_string(const T &value,
+ const Patterns::PatternBase &p = *Convert<T>::to_pattern())
+ {
+ namespace B = magic_enum::bitwise_operators;
+ const auto values = magic_enum::enum_values<T>();
+ std::vector<std::string> names;
+ for (const auto &v : values)
+ if (B::operator&(value, v) == v)
+ names.push_back(std::string(magic_enum::enum_name(v)));
+ return Patterns::Tools::Convert<decltype(names)>::to_string(names, p);
+ }
+
+ static T
+ to_value(const std::string &s,
+ const dealii::Patterns::PatternBase &p = *to_pattern())
+ {
+ namespace B = magic_enum::bitwise_operators;
+ // Make sure we have a valid enum value, or empty value
+ AssertThrow(p.match(s), ExcNoMatch(s, p.description()));
+ T value = T();
+ std::vector<std::string> value_names;
+ value_names =
+ Patterns::Tools::Convert<decltype(value_names)>::to_value(s, p);
+ for (const auto &name : value_names)
+ {
+ auto v = magic_enum::enum_cast<T>(name);
+ if (v.has_value())
+ value = B::operator|(value, v.value());
+ }
+ return value;
+ }
+ };
+#endif
+
// Utility function with default Pattern
template <typename T>
std::string
--- /dev/null
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2010 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+// Check Patterns::Tools::Convert for enum types
+
+#include <deal.II/fe/fe_values.h>
+
+#include <memory>
+
+#include "../tests.h"
+
+int
+main()
+{
+ initlog();
+
+ UpdateFlags flags = update_values;
+
+ deallog << Patterns::Tools::Convert<UpdateFlags>::to_string(flags)
+ << std::endl;
+
+ flags = Patterns::Tools::Convert<UpdateFlags>::to_value(
+ "update_values|update_gradients");
+
+ deallog << Patterns::Tools::Convert<UpdateFlags>::to_string(flags)
+ << std::endl;
+}
--- /dev/null
+
+DEAL::update_default| update_values
+DEAL::update_default| update_values| update_gradients