From: Johannes Heinz Date: Mon, 8 Jan 2024 20:15:48 +0000 (+0100) Subject: rename template parameters X-Git-Tag: relicensing~187^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=b71ee361544ccf72cfec0cb2fc83fcf9f73ac6db;p=dealii.git rename template parameters --- diff --git a/include/deal.II/cgal/intersections.h b/include/deal.II/cgal/intersections.h index 1fa1236682..f9ef1757e7 100644 --- a/include/deal.II/cgal/intersections.h +++ b/include/deal.II/cgal/intersections.h @@ -50,14 +50,14 @@ namespace CGALWrappers * @param tol Threshold to decide whether or not a simplex is included. * @return Vector of arrays, where each array identify a simplex by its vertices. */ - template - std::vector, dim1 + 1>> + template + std::vector, structdim1 + 1>> compute_intersection_of_cells( - const typename Triangulation::cell_iterator &cell0, - const typename Triangulation::cell_iterator &cell1, - const Mapping &mapping0, - const Mapping &mapping1, - const double tol = 1e-9); + const typename Triangulation::cell_iterator &cell0, + const typename Triangulation::cell_iterator &cell1, + const Mapping &mapping0, + const Mapping &mapping1, + const double tol = 1e-9); /** @@ -68,8 +68,8 @@ namespace CGALWrappers * * @note The vertices have to be given in CGAL order. */ - template - std::vector, dim1 + 1>> + template + std::vector, structdim1 + 1>> compute_intersection_of_cells( const ArrayView> &vertices0, const ArrayView> &vertices1, diff --git a/source/cgal/intersections.cc b/source/cgal/intersections.cc index 48bfe32fb5..a75999eaf1 100644 --- a/source/cgal/intersections.cc +++ b/source/cgal/intersections.cc @@ -749,8 +749,8 @@ namespace CGALWrappers } // namespace internal - template - std::vector, dim1 + 1>> + template + std::vector, structdim1 + 1>> compute_intersection_of_cells( const ArrayView> &vertices0, const ArrayView> &vertices1, @@ -764,7 +764,7 @@ namespace CGALWrappers ExcMessage( "The intersection cannot be computed as at least one of the two cells has no vertices.")); - if constexpr (dim0 == 2 && dim1 == 2 && spacedim == 2) + if constexpr (structdim0 == 2 && structdim1 == 2 && spacedim == 2) { if (n_vertices0 == 4 && n_vertices1 == 4) { @@ -773,7 +773,7 @@ namespace CGALWrappers tol); } } - else if constexpr (dim0 == 2 && dim1 == 1 && spacedim == 2) + else if constexpr (structdim0 == 2 && structdim1 == 1 && spacedim == 2) { if (n_vertices0 == 4 && n_vertices1 == 2) { @@ -782,7 +782,7 @@ namespace CGALWrappers tol); } } - else if constexpr (dim0 == 3 && dim1 == 1 && spacedim == 3) + else if constexpr (structdim0 == 3 && structdim1 == 1 && spacedim == 3) { if (n_vertices0 == 8 && n_vertices1 == 2) { @@ -791,7 +791,7 @@ namespace CGALWrappers tol); } } - else if constexpr (dim0 == 3 && dim1 == 2 && spacedim == 3) + else if constexpr (structdim0 == 3 && structdim1 == 2 && spacedim == 3) { if (n_vertices0 == 8 && n_vertices1 == 4) { @@ -800,7 +800,7 @@ namespace CGALWrappers tol); } } - else if constexpr (dim0 == 3 && dim1 == 3 && spacedim == 3) + else if constexpr (structdim0 == 3 && structdim1 == 3 && spacedim == 3) { if (n_vertices0 == 8 && n_vertices1 == 8) { @@ -819,20 +819,20 @@ namespace CGALWrappers } - template - std::vector, dim1 + 1>> + template + std::vector, structdim1 + 1>> compute_intersection_of_cells( - const typename Triangulation::cell_iterator &cell0, - const typename Triangulation::cell_iterator &cell1, - const Mapping &mapping0, - const Mapping &mapping1, - const double tol) + const typename Triangulation::cell_iterator &cell0, + const typename Triangulation::cell_iterator &cell1, + const Mapping &mapping0, + const Mapping &mapping1, + const double tol) { Assert(mapping0.get_vertices(cell0).size() == - ReferenceCells::get_hypercube().n_vertices(), + ReferenceCells::get_hypercube().n_vertices(), ExcNotImplemented()); Assert(mapping1.get_vertices(cell1).size() == - ReferenceCells::get_hypercube().n_vertices(), + ReferenceCells::get_hypercube().n_vertices(), ExcNotImplemented()); const auto &vertices0 = @@ -840,9 +840,8 @@ namespace CGALWrappers const auto &vertices1 = CGALWrappers::get_vertices_in_cgal_order(cell1, mapping1); - return compute_intersection_of_cells(vertices0, - vertices1, - tol); + return compute_intersection_of_cells( + vertices0, vertices1, tol); } # include "intersections.inst" @@ -855,12 +854,12 @@ DEAL_II_NAMESPACE_CLOSE DEAL_II_NAMESPACE_OPEN -template -std::vector, dim1 + 1>> +std::vector, structdim1 + 1>> compute_intersection_of_cells( const std::array, n_components0> &vertices0, const std::array, n_components1> &vertices1, @@ -872,14 +871,14 @@ compute_intersection_of_cells( AssertThrow(false, ExcNeedsCGAL()); } -template -std::vector, dim1 + 1>> +template +std::vector, structdim1 + 1>> compute_intersection_of_cells( - const typename Triangulation::cell_iterator &cell0, - const typename Triangulation::cell_iterator &cell1, - const Mapping &mapping0, - const Mapping &mapping1, - const double tol) + const typename Triangulation::cell_iterator &cell0, + const typename Triangulation::cell_iterator &cell1, + const Mapping &mapping0, + const Mapping &mapping1, + const double tol) { (void)cell0; (void)cell1; diff --git a/source/cgal/intersections.inst.in b/source/cgal/intersections.inst.in index 47b1327a8b..98de036407 100644 --- a/source/cgal/intersections.inst.in +++ b/source/cgal/intersections.inst.in @@ -15,20 +15,21 @@ -for (dim0 : DIMENSIONS; dim1 : DIMENSIONS; spacedim : SPACE_DIMENSIONS) +for (structdim0 : DIMENSIONS; structdim1 : DIMENSIONS; + spacedim : SPACE_DIMENSIONS) { -#if dim0 <= spacedim && dim0 >= dim1 +#if structdim0 <= spacedim && structdim0 >= structdim1 - template std::vector, dim1 + 1>> - compute_intersection_of_cells( - const typename Triangulation::cell_iterator &cell0, - const typename Triangulation::cell_iterator &cell1, - const Mapping &mapping0, - const Mapping &mapping1, - const double tol); + template std::vector, structdim1 + 1>> + compute_intersection_of_cells( + const typename Triangulation::cell_iterator &cell0, + const typename Triangulation::cell_iterator &cell1, + const Mapping &mapping0, + const Mapping &mapping1, + const double tol); - template std::vector, dim1 + 1>> - compute_intersection_of_cells( + template std::vector, structdim1 + 1>> + compute_intersection_of_cells( const ArrayView> &vertices0, const ArrayView> &vertices1, const double tol);