]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add the scripts that convert files to modules.
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 23 May 2025 15:11:02 +0000 (09:11 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Tue, 27 May 2025 14:56:18 +0000 (08:56 -0600)
contrib/utilities/convert_header_file_to_interface_module_unit.py [new file with mode: 0644]
contrib/utilities/convert_source_file_to_implementation_module_unit.py [new file with mode: 0644]
contrib/utilities/convert_to_module_units_common.py [new file with mode: 0644]

diff --git a/contrib/utilities/convert_header_file_to_interface_module_unit.py b/contrib/utilities/convert_header_file_to_interface_module_unit.py
new file mode 100644 (file)
index 0000000..f20286f
--- /dev/null
@@ -0,0 +1,136 @@
+#!/usr/bin/python3
+
+## ------------------------------------------------------------------------
+##
+## SPDX-License-Identifier: LGPL-2.1-or-later
+## Copyright (C) 2025 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.
+##
+## ------------------------------------------------------------------------
+
+# Given a C++ header file as input, read through it and output an
+# equivalent interface module unit that exports a module partition
+# with the declarations that are part of the header file, and that
+# uses 'import' statements in place of existing '#include' directives.
+#
+# Call this script via
+#   python3 contrib/utilities/convert_header_file_to_interface_module_unit.py in.h out.ccm
+
+
+import sys
+import re
+
+import convert_to_module_units_common as header_to_partition_maps
+
+
+header_file                = sys.argv[1]
+interface_module_unit_file = sys.argv[2]
+
+assert header_file, "No input file name given."
+assert interface_module_unit_file, "No output file name given."
+
+m = re.search(r".*/include/deal.II/(.*).h", header_file)
+interface_module_partition_name = "interface_partition_" + m.group(1).replace("/", "_")
+
+
+# Read the entire header file into memory:
+input = open(header_file, "r")
+lines = input.readlines()
+input.close()
+
+
+
+# Regular expressions that match the preprocessor defines that we use
+# to open and close the deal.II namespace. We allow for whitespace
+# before and after, but nothing else in these regular
+# expressions. That means that the lines in question cannot contain
+# anything else -- which we will use in a few places to ensure that
+# these statements are not expanded.
+match_dealii_open = re.compile(r"^ *DEAL_II_NAMESPACE_OPEN *$")
+match_dealii_close = re.compile(r"^ *DEAL_II_NAMESPACE_CLOSE *$")
+
+# Then scan through it and collect the deal.II headers that are
+# '#include'd, excluding "config.h" and "exception_macros.h", both of
+# which declare macros and so much be #included rather than imported:
+dealii_include_list = header_to_partition_maps.get_dealii_includes(lines)
+
+# Determine which external projects that we wrap with module
+# partitions are being utilized by the current file:
+used_external_projects = header_to_partition_maps.get_used_external_projects(lines)
+
+
+# Now open the output file and start writing into it. Start the file
+# with the 'module;' text that indicates that this is going to be a
+# module file. Then also include the file that defines all deal.II
+# macros (including everything that's in config.h) that ensures that
+# we can use preprocessor defines describing the configuration of
+# deal.II. We need to explicitly include it here because we comment
+# out all other #includes and so can't get at config.h transitively.
+output = open(interface_module_unit_file, "w")
+output.write("module;\n")
+output.write("#include <deal.II/macros.h>\n\n")
+
+# Then go through the lines of the input file again:
+for line in lines :
+
+    # If this was an '#include' directive for a deal.II header, or an
+    # external project that we wrap, then just ignore/remove this
+    # line.
+    #
+    # In source files, it is in principle acceptable to have #includes
+    # because source files do not have exported partitions and so
+    # processing #includes does not contribute to the bloat that
+    # happens when importing interface partitions that each contain
+    # dozens or hundreds of (transitive) includes. As a consequence,
+    # we do allow #includes in specific circumstances, when they are
+    # specifically marked with a comment. But we do not allow this
+    # in interface units (header files) and so error out if someone
+    # tries this.
+    if (header_to_partition_maps.match_dealii_includes.match(line) or
+        header_to_partition_maps.matches_external_header_include(line)) :
+        if "Do not convert for module purposes" in line :
+            raise "It is not allowed to escape the conversion of include statements in header files."
+        else :
+            pass
+
+    # If this line contained the text DEAL_II_NAMESPACE_OPEN, then
+    # prefix this text by the interface module unit start, the import
+    # statements that correspond to the previously '#include'd header
+    # files, and an 'export {' statement.
+    #
+    # Because import statements do not transitively import what the
+    # imported module itself imported, we would have to annotate every
+    # source file with a complete list of all module partitions (or,
+    # better, corresponding header files) it needs something
+    # from. This is not realistic. Rather, we emulate the semantics of
+    # header files by not only importing partitions, but re-exporting
+    # them as well via the 'export import :interface_partition;'
+    # statement.
+    elif match_dealii_open.match(line) :
+        output.write("\nexport module dealii : " + interface_module_partition_name + ";\n\n")
+
+        for external_project in used_external_projects :
+            output.write("import :interface_partition_" + external_project + ";\n")
+
+        for inc in dealii_include_list :
+            module_partition = "interface_partition_" + inc.replace("/", "_")
+            output.write("import :" + module_partition + ";\n")
+
+        output.write("\nexport\n{\n\n")
+        output.write(line)  # Copy the previous line
+
+    # If this line contained the text DEAL_II_NAMESPACE_CLODE, then
+    # close the previous 'export {' statement:
+    elif match_dealii_close.match(line) :
+        output.write(line)  # Copy the previous line
+        output.write("\n} // close previous 'export' statement\n\n")
+
+    # Otherwise just copy the previous line:
+    else :
+        output.write(line)
diff --git a/contrib/utilities/convert_source_file_to_implementation_module_unit.py b/contrib/utilities/convert_source_file_to_implementation_module_unit.py
new file mode 100644 (file)
index 0000000..8e95a89
--- /dev/null
@@ -0,0 +1,112 @@
+#!/usr/bin/python3
+
+## ------------------------------------------------------------------------
+##
+## SPDX-License-Identifier: LGPL-2.1-or-later
+## Copyright (C) 2025 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.
+##
+## ------------------------------------------------------------------------
+
+# Given a C++ source file (a .cc file) as input, read through it and
+# output an equivalent implementation module unit that provides a
+# (non-exported) module partition with the stuff that is part of the
+# header file, and that uses 'import' statements in place of existing
+# '#include' directives.
+#
+# Call this script via
+#   python3 contrib/utilities/convert_source_file_to_implementation_module_unit.py in.cc out.cc
+
+
+import sys
+import re
+
+import convert_to_module_units_common as header_to_partition_maps
+
+
+source_file                = sys.argv[1]
+implementation_module_unit_file = sys.argv[2]
+
+assert source_file, "No input file name given."
+assert implementation_module_unit_file, "No output file name given."
+
+m = re.search(r".*source/(.*).cc", source_file)
+implementation_module_partition_name = "implementation_partition_" + m.group(1).replace("/", "_")
+
+
+# Read the entire source file into memory:
+input = open(source_file, "r")
+lines = input.readlines()
+input.close()
+
+
+match_dealii_open = re.compile(r" *DEAL_II_NAMESPACE_OPEN *")
+
+# Then scan through it and collect the deal.II headers that are
+# '#include'd, excluding "config.h" and "exception_macros.h", both of
+# which declare macros and so much be #included rather than imported:
+dealii_include_list = header_to_partition_maps.get_dealii_includes(lines)
+
+# Determine which external projects that we wrap with module
+# partitions are being utilized by the current file:
+used_external_projects = header_to_partition_maps.get_used_external_projects(lines)
+
+
+# Now open the output file and start writing into it. Start the file
+# with the 'module;' text that indicates that this is going to be a
+# module file. Then also include the file that defines all deal.II
+# macros (including everything that's in config.h) that ensures that
+# we can use preprocessor defines describing the configuration of
+# deal.II. We need to explicitly include it here because we comment
+# out all other #includes and so can't get at config.h transitively.
+output = open(implementation_module_unit_file, "w")
+output.write("module;\n")
+output.write("#include <deal.II/macros.h>\n\n")
+
+
+# Then go through the lines of the input file again:
+for line in lines :
+
+    # If this was an '#include' directive for a deal.II header, or an
+    # external project that we wrap, then just ignore/remove this
+    # line.
+    #
+    # In source files, it is in principle acceptable to have #includes
+    # because source files do not have exported partitions and so
+    # processing #includes does not contribute to the bloat that
+    # happens when importing interface partitions that each contain
+    # dozens or hundreds of (transitive) includes. As a consequence,
+    # we do allow #includes in specific circumstances, when they are
+    # specifically marked with a comment.
+    if (header_to_partition_maps.match_dealii_includes.match(line) or
+        header_to_partition_maps.matches_external_header_include(line)) :
+        if "Do not convert for module purposes" in line :
+            output.write(line)
+        else :
+            pass
+
+    # If this line contained the text DEAL_II_NAMESPACE_OPEN, then
+    # prefix this text by the implementation module unit start, the import
+    # statements that correspond to the previously '#include'd header
+    # files, but importantly no 'export {' statement:
+    elif match_dealii_open.match(line) :
+        output.write("module dealii : " + implementation_module_partition_name + ";\n\n")
+
+        for external_project in used_external_projects :
+            output.write("import :interface_partition_" + external_project + ";\n")
+
+        for inc in dealii_include_list :
+            module_partition = "interface_partition_" + inc.replace("/", "_")
+            output.write("import :" + module_partition + ";\n")
+        output.write("\n")
+        output.write(line)  # Copy the previous line
+
+    # Otherwise just copy the previous line:
+    else :
+        output.write(line)
diff --git a/contrib/utilities/convert_to_module_units_common.py b/contrib/utilities/convert_to_module_units_common.py
new file mode 100644 (file)
index 0000000..d385640
--- /dev/null
@@ -0,0 +1,341 @@
+## ------------------------------------------------------------------------
+##
+## SPDX-License-Identifier: LGPL-2.1-or-later
+## Copyright (C) 2025 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.
+##
+## ------------------------------------------------------------------------
+
+
+# This file defines some common infrastructure for the two scripts
+# that convert header and source files into module partition
+# units. Specifically, it contains lists of header files for each
+# external dependency of deal.II, and compiles them into regular
+# expressions that can be matched against #include statements.
+
+
+import re
+
+# Define lists of (potentially wildcarded) names of header files for
+# deal.II and all other external dependencies we use.. The list of C++
+# headers is from
+# https://stackoverflow.com/questions/2027991/list-of-standard-header-files-in-c-and-c
+dealii_includes = [ "deal.II/(.*).h" ]
+
+adolc_includes  = [ "adolc/.*" ]
+
+boost_includes  = [ "boost/.*" ]
+
+cgal_includes  = [ "CGAL/.*" ]
+
+hdf5_includes  = [ "hdf5.h" ]
+
+kokkos_includes = [ "Kokkos.*" ]
+
+metis_includes = [ "metis.h" ]
+
+mpi_includes = [ "mpi.h" ]
+
+mumps_includes = [ "mumps_.*.h" ]
+
+muparser_includes = [ "muParser.h" ]
+
+opencascade_includes = [ "Adaptor3d_Curve.hxx",
+                         "Adaptor3d_HCurve.hxx",
+                         "BRep.*.hxx",
+                         "Geom.*.hxx",
+                         "Handle.*.hxx",
+                         "IGESControl.*.hxx",
+                         "IFSelect_ReturnStatus.hxx",
+                         "IntCurve.*.hxx",
+                         "Poly_Triangulation.hxx",
+                         "ShapeAnalysis.*.hxx",
+                         "Standard_Transient.hxx",
+                         "STEPControl.*.hxx",
+                         "StlAPI.*.hxx",
+                         "TCol.*.hxx",
+                         "Top.*.hxx",
+                         "TopoDS.*.hxx",
+                         "gp_.*.hxx" ]
+
+p4est_includes  = [ "p4est.*",
+                    "p8est.*",
+                    "sc_containers.h" ]
+
+slepc_includes = [ "slepc.*.h" ]
+
+petsc_includes  = [ "petscconf.h",
+                    "petscdm.h",
+                    "petscerror.h",
+                    "petscis.h",
+                    "petscistypes.h",
+                    "petscksp.h",
+                    "petscmat.h",
+                    "petscpc.h",
+                    "petscsf.h",
+                    "petscsnes.h",
+                    "petscsys.h",
+                    "petscts.h",
+                    "petscvec.h",
+                    "petsc/.*.h" ]
+
+std_includes    = [ "assert.h",
+                    "limits.h",
+                    "signal.h",
+                    "stdlib.h",
+                    "ctype.h",
+                    "locale.h",
+                    "stdarg.h",
+                    "string.h",
+                    "errno.h",
+                    "math.h",
+                    "stddef.h",
+                    "time.h",
+                    "float.h",
+                    "setjmp.h",
+                    "stdio.h",
+                    "algorithm",
+                    "format",
+                    "new",
+                    "stdexcept",
+                    "any",
+                    "forward_list",
+                    "numbers",
+                    "stdfloat",
+                    "array",
+                    "fstream",
+                    "numeric",
+                    "stop_token",
+                    "atomic",
+                    "functional",
+                    "optional",
+                    "streambuf",
+                    "barrier",
+                    "future",
+                    "ostream",
+                    "string",
+                    "bit",
+                    "generator",
+                    "print",
+                    "string_view",
+                    "bitset",
+                    "initializer_list",
+#                    "queue",   # excluded for now due to https://github.com/llvm/llvm-project/issues/138558
+                    "strstream",
+                    "charconv",
+                    "iomanip",
+#                    "random",   # excluded for now due to https://github.com/llvm/llvm-project/issues/138558
+                    "syncstream",
+                    "chrono",
+                    "ios",
+                    "ranges",
+                    "system_error",
+                    "codecvt",
+                    "iosfwd",
+                    "ratio",
+                    "thread",
+                    "compare",
+                    "iostream",
+                    "regex",
+                    "tuple",
+                    "complex",
+                    "istream",
+                    "scoped_allocator",
+                    "type_traits",
+                    "concepts",
+                    "iterator",
+                    "semaphore",
+                    "typeindex",
+                    "condition_variable",
+                    "latch",
+                    "set",
+                    "typeinfo",
+                    "coroutine",
+                    "limits",
+                    "shared_mutex",
+                    "unordered_map",
+                    "deque",
+                    "list",
+                    "source_location",
+                    "unordered_set",
+                    "exception",
+                    "locale",
+                    "span",
+                    "utility",
+                    "execution",
+                    "map",
+                    "spanstream",
+                    "valarray",
+                    "expected",
+                    "mdspan",
+                    "sstream",
+                    "variant",
+                    "filesystem",
+                    "memory",
+#                    "stack",   # excluded for now due to https://github.com/llvm/llvm-project/issues/138558
+                    "vector",
+                    "flat_map",
+                    "memory_resource",
+                    "stacktrace",
+                    "version",
+                    "flat_set",
+                    "mutex",
+                    "cassert",
+                    "cfenv",
+                    "climits",
+                    "csetjmp",
+                    "cstddef",
+                    "cstdlib",
+                    "cuchar",
+                    "cctype",
+                    "cfloat",
+                    "clocale",
+                    "csignal",
+                    "cstdint",
+                    "cstring",
+                    "cwchar",
+                    "cerrno",
+                    "cinttypes",
+                    "cmath",
+                    "cstdarg",
+                    "cstdio",
+                    "ctime",
+                    "cwctype" ]
+
+sundials_includes = [ "arkode/arkode.*",
+                      "ida/ida.h",
+                      "idas/idas.h",
+                      "kinsol/kinsol.*",
+                      "nvector/nvector_parallel.h",
+                      "nvector/nvector_serial.h",
+                      "sundials/sundials_.*",
+                      "sunlinsol/sunlinsol_.*",
+                      "sunmatrix/sunmatrix_.*",
+                      "sunnonlinsol/sunnonlinsol_.*" ]
+
+taskflow_includes = [ "taskflow/.*" ]
+
+tbb_includes = [ "tbb/.*" ]
+
+trilinos_includes = [ "Amesos.*",
+                      "AztecOO.h",
+                      "Belos.*",
+                      "Epetra.*",
+                      "Ifpack2.*",
+                      "NOX.*",
+                      "Rol.*",
+                      "ROL.*",
+                      "Sacado.*",
+                      "Teuchos.*",
+                      "Tpetra.*",
+                      "zoltan_cpp.h",
+                      "exodusII.h" ]
+
+umfpack_includes = [ "umfpack.h" ]
+
+vtk_includes = [ "vtk.*" ]
+
+zlib_includes = [ "zlib.h" ]
+
+# Define regexs to match header #includes for each of the categories above
+match_dealii_includes      = re.compile("# *include *<" + "|".join(dealii_includes) + ">")
+
+match_adolc_includes       = re.compile("# *include *<(" + "|".join(adolc_includes)      + ")>")
+match_boost_includes       = re.compile("# *include *<(" + "|".join(boost_includes)      + ")>")
+match_cgal_includes        = re.compile("# *include *<(" + "|".join(cgal_includes)       + ")>")
+match_hdf5_includes        = re.compile("# *include *<(" + "|".join(hdf5_includes)       + ")>")
+match_kokkos_includes      = re.compile("# *include *<(" + "|".join(kokkos_includes)     + ")>")
+match_metis_includes       = re.compile("# *include *<(" + "|".join(metis_includes)      + ")>")
+match_mpi_includes         = re.compile("# *include *<(" + "|".join(mpi_includes)        + ")>")
+match_mumps_includes       = re.compile("# *include *<(" + "|".join(mumps_includes)        + ")>")
+match_muparser_includes    = re.compile("# *include *<(" + "|".join(muparser_includes)        + ")>")
+match_opencascade_includes = re.compile("# *include *<(" + "|".join(opencascade_includes)+ ")>")
+match_p4est_includes       = re.compile("# *include *<(" + "|".join(p4est_includes)      + ")>")
+match_petsc_includes       = re.compile("# *include *<(" + "|".join(petsc_includes)      + ")>")
+match_slepc_includes       = re.compile("# *include *<(" + "|".join(slepc_includes)      + ")>")
+match_std_includes         = re.compile("# *include *<(" + "|".join(std_includes)        + ")>")
+match_sundials_includes    = re.compile("# *include *<(" + "|".join(sundials_includes)   + ")>")
+match_taskflow_includes    = re.compile("# *include *<(" + "|".join(taskflow_includes)        + ")>")
+match_tbb_includes         = re.compile("# *include *<(" + "|".join(tbb_includes)        + ")>")
+match_trilinos_includes    = re.compile("# *include *<(" + "|".join(trilinos_includes)   + ")>")
+match_umfpack_includes     = re.compile("# *include *<(" + "|".join(umfpack_includes)    + ")>")
+match_vtk_includes         = re.compile("# *include *<(" + "|".join(vtk_includes)        + ")>")
+match_zlib_includes        = re.compile("# *include *<(" + "|".join(zlib_includes)       + ")>")
+
+external_package_headers_regex_map = { "adolc"       : match_adolc_includes,
+                                       "boost"       : match_boost_includes,
+                                       "cgal"        : match_cgal_includes,
+                                       "hdf5"        : match_hdf5_includes,
+                                       "kokkos"      : match_kokkos_includes,
+                                       "metis"       : match_metis_includes,
+                                       "mpi"         : match_mpi_includes,
+                                       "mumps"       : match_mumps_includes,
+                                       "muparser"    : match_muparser_includes,
+                                       "opencascade" : match_opencascade_includes,
+                                       "p4est"       : match_p4est_includes,
+                                       "petsc"       : match_petsc_includes,
+                                       "slepc"       : match_slepc_includes,
+                                       "std"         : match_std_includes,
+                                       "sundials"    : match_sundials_includes,
+                                       "taskflow"    : match_taskflow_includes,
+                                       "tbb"         : match_tbb_includes,
+                                       "trilinos"    : match_trilinos_includes,
+                                       "umfpack"     : match_umfpack_includes,
+                                       "vtk"         : match_vtk_includes,
+                                       "zlib"        : match_zlib_includes }
+
+
+# A function that given a list of strings (specifically: the
+# individual lines of a file) returns the list of deal.II header files
+# that are being #included, excluding "config.h" and "exception_macros.h".
+def get_dealii_includes(lines) :
+    dealii_include_list = []
+    for line in lines :
+        m = match_dealii_includes.match(line)
+        if m :
+            if (m.group(1) != "base/config") and (m.group(1) != "base/exception_macros") :
+                dealii_include_list.append(m.group(1))
+    return dealii_include_list
+
+
+
+# Return whether a line in a header or source file matches an #include
+# statement for any of the recognized external projects that we wrap
+# in module partitions:
+def matches_external_header_include(line) :
+    for package,regex in external_package_headers_regex_map.items() :
+        if regex.match(line) :
+            return True
+    return False
+
+
+# A function that given a list of strings (specifically: the
+# individual lines of a file) returns which of the external projects
+# we wrap in module partitions are actually being used by this one
+# file (as indicated by this file #including any of the corresponding
+# header files of that external project).
+def get_used_external_projects(lines) :
+    used_external_projects = set()
+    my_external_package_headers_regex_map = external_package_headers_regex_map.copy()
+    
+    for line in lines :
+        for package, regex in my_external_package_headers_regex_map.items() :
+            if regex.match(line) :
+                # The current line matches an external package's
+                # headers. Once we have determined that a file uses an
+                # external package, we no longer need to test for that and
+                # can remove that entry from the list of regexes. Note that
+                # for this to work, we need to use .copy() above when creating
+                # the my_... dictionary, as we would otherwise simply be
+                # deleting the entry in the original dictionary.
+                used_external_projects.add(package)
+                del my_external_package_headers_regex_map[package]
+                break
+
+    return used_external_projects

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.