]> https://gitweb.dealii.org/ - dealii.git/commitdiff
fix: building problem with MPI on Windows system 13593/head
authorDaniel Sun <sunnyboypdan@126.com>
Tue, 5 Apr 2022 15:53:29 +0000 (23:53 +0800)
committerDaniel Sun <sunnyboypdan@126.com>
Tue, 26 Apr 2022 08:05:51 +0000 (16:05 +0800)
doc/news/changes/minor/20220426DanielSun [new file with mode: 0644]
include/deal.II/fe/fe_tools_interpolate.templates.h
include/deal.II/grid/grid_tools.h
include/deal.II/grid/grid_tools_cache.h
include/deal.II/numerics/vector_tools_constraints.h
include/deal.II/particles/generators.h

diff --git a/doc/news/changes/minor/20220426DanielSun b/doc/news/changes/minor/20220426DanielSun
new file mode 100644 (file)
index 0000000..344ecae
--- /dev/null
@@ -0,0 +1,5 @@
+Fixed: The building process now works well with MPI option on Windows system.
+If configured properly, we can explore the capability of distributed computing
+on Windows system within Deal.II framework.
+<br>
+(Daniel Sun, 2022/04/26)
index cac9b99d7496e83123c7b742710dfa6eb72ef7fd..fbdf6ff1b6f35986fda2b78ccdafb80c65484f04 100644 (file)
@@ -738,7 +738,7 @@ namespace FETools
       interpolation_difference(dof1, u1, dof2.get_fe(), u1_difference);
     else
       {
-        internal::interpolation_difference(
+        internal::interpolation_difference<dim, InVector, OutVector, spacedim>(
           dof1, constraints1, u1, dof2, constraints2, u1_difference);
       }
   }
index 8bf837ebd7b4d59eae68bb654c89b6fcc1219dae..2ac4cf141b44a789781ba1dce93a56d7733c93bc 100644 (file)
@@ -176,7 +176,12 @@ namespace GridTools
   volume(const Triangulation<dim, spacedim> &tria,
          const Mapping<dim, spacedim> &      mapping =
            (ReferenceCells::get_hypercube<dim>()
-              .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+              .template get_default_linear_mapping<dim, spacedim>()
+#else
+              .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+              ));
 
   /**
    * Return an approximation of the diameter of the smallest active cell of a
@@ -194,7 +199,12 @@ namespace GridTools
     const Triangulation<dim, spacedim> &triangulation,
     const Mapping<dim, spacedim> &      mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Return an approximation of the diameter of the largest active cell of a
@@ -212,7 +222,12 @@ namespace GridTools
     const Triangulation<dim, spacedim> &triangulation,
     const Mapping<dim, spacedim> &      mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Given a list of vertices (typically obtained using
@@ -1243,7 +1258,12 @@ namespace GridTools
     const Triangulation<dim, spacedim> &container,
     const Mapping<dim, spacedim> &      mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Find and return the index of the closest vertex to a given point in the
@@ -2069,7 +2089,12 @@ namespace GridTools
     const Point<spacedim> &                                            position,
     const Mapping<dim, spacedim> &                                     mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Compute a globally unique index for each vertex and hanging node
@@ -4358,10 +4383,13 @@ namespace GridTools
       for (const auto ghost_owner : ghost_owners)
         neighbor_cell_list[ghost_owner] = {};
 
-      process_cells([&](const auto &cell, const auto key) {
+      process_cells([&](const auto &cell, const auto key) -> void {
         if (cell_filter(cell))
-          neighbor_cell_list[key].emplace_back(
-            cell->id().template to_binary<spacedim>());
+          {
+            constexpr int spacedim = MeshType::space_dimension;
+            neighbor_cell_list[key].emplace_back(
+              cell->id().template to_binary<spacedim>());
+          }
       });
 
       Assert(ghost_owners.size() == neighbor_cell_list.size(),
index 62db3395d1db8cf741f2e9db22a5a34e5dbc6c68..f188fd791cfd7b2fa038894a7b618560f9e6be66 100644 (file)
@@ -77,7 +77,12 @@ namespace GridTools
     Cache(const Triangulation<dim, spacedim> &tria,
           const Mapping<dim, spacedim> &      mapping =
             (ReferenceCells::get_hypercube<dim>()
-               .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+               .template get_default_linear_mapping<dim, spacedim>()
+#else
+               .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+               ));
 
     /**
      * Destructor.
index ff7b6a7a3c2dada2be7c5c10ff4e9dd5f596eed1..cb29094c8b5f9bfb8512d97bfdeda5d937dc8707 100644 (file)
@@ -281,7 +281,12 @@ namespace VectorTools
     AffineConstraints<double> &   constraints,
     const Mapping<dim, spacedim> &mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * This function does the same as the
@@ -304,7 +309,12 @@ namespace VectorTools
     AffineConstraints<double> &         constraints,
     const Mapping<dim, spacedim> &      mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Compute the constraints that correspond to boundary conditions of the
@@ -333,7 +343,12 @@ namespace VectorTools
     AffineConstraints<double> &   constraints,
     const Mapping<dim, spacedim> &mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   /**
    * Same as above for homogeneous tangential-flux constraints.
@@ -352,7 +367,12 @@ namespace VectorTools
     AffineConstraints<double> &         constraints,
     const Mapping<dim, spacedim> &      mapping =
       (ReferenceCells::get_hypercube<dim>()
-         .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+         .template get_default_linear_mapping<dim, spacedim>()
+#else
+         .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+         ));
 
   //@}
 } // namespace VectorTools
index 5059bf415d7ec28fa853de9c3cb74d935c5f7825..20dc792231d4ee654e2158d02d8fa08b108e4db9 100644 (file)
@@ -70,7 +70,12 @@ namespace Particles
       ParticleHandler<dim, spacedim> &    particle_handler,
       const Mapping<dim, spacedim> &      mapping =
         (ReferenceCells::get_hypercube<dim>()
-           .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+           .template get_default_linear_mapping<dim, spacedim>()
+#else
+           .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+           ));
 
     /**
      * A function that generates one particle at a random location in cell @p cell and with
@@ -109,7 +114,12 @@ namespace Particles
       std::mt19937 &                random_number_generator,
       const Mapping<dim, spacedim> &mapping =
         (ReferenceCells::get_hypercube<dim>()
-           .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+           .template get_default_linear_mapping<dim, spacedim>()
+#else
+           .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+           ));
 
     /**
      * A function that generates one particle at a random location in cell @p cell and with
@@ -126,7 +136,12 @@ namespace Particles
       ParticleHandler<dim, spacedim> &particle_handler,
       const Mapping<dim, spacedim> &  mapping =
         (ReferenceCells::get_hypercube<dim>()
-           .template get_default_linear_mapping<dim, spacedim>()));
+#ifndef _MSC_VER
+           .template get_default_linear_mapping<dim, spacedim>()
+#else
+           .ReferenceCell::get_default_linear_mapping<dim, spacedim>()
+#endif
+           ));
 
     /**
      * A function that generates particles randomly in the domain with a
@@ -182,7 +197,11 @@ namespace Particles
       ParticleHandler<dim, spacedim> &    particle_handler,
       const Mapping<dim, spacedim> &      mapping =
         (ReferenceCells::get_hypercube<dim>()
+#ifndef _MSC_VER
            .template get_default_linear_mapping<dim, spacedim>()),
+#else
+           .ReferenceCell::get_default_linear_mapping<dim, spacedim>()),
+#endif
       const unsigned int random_number_seed = 5432);
 
 
@@ -229,7 +248,11 @@ namespace Particles
       ParticleHandler<dim, spacedim> &particle_handler,
       const Mapping<dim, spacedim> &  mapping =
         (ReferenceCells::get_hypercube<dim>()
+#ifndef _MSC_VER
            .template get_default_linear_mapping<dim, spacedim>()),
+#else
+           .ReferenceCell::get_default_linear_mapping<dim, spacedim>()),
+#endif
       const ComponentMask &                   components = ComponentMask(),
       const std::vector<std::vector<double>> &properties = {});
 

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.