MFEM  v4.6.0
Finite element discretization library
nodal-transfer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details
11 //
12 // -------------------------------------------------------------------------
13 // Nodal Transfer Miniapp: Map ParGridFunction to Different MPI Partitioning
14 // -------------------------------------------------------------------------
15 //
16 // The Nodal Transfer Miniapp maps partitioned parallel grid function to a
17 // parallel grid function partitioned on a different number of processes. The
18 // miniapp has two regimes: 1) Generates partitioned parallel grid function
19 // and saves it to a set of files; 2) Reads the partitioned grid function and
20 // maps it to the current partition. The map assumes that the position of the
21 // nodal DOFs does not change between the original grid function and the target
22 // grid function. The transfer does not perform any interpolation. It just
23 // copies the nodal values between the two grid functions.
24 //
25 // Generate second order mesh on 4 processes
26 // mpirun -np 4 ./nodal-transfer -rs 2 -rp 1 -gd 1 -o 2
27 // Read the generated data and map it to a grid function defined on two processes
28 // mpirun -np 2 ./nodal-transfer -rs 2 -rp 0 -gd 0 -snp 4 -o 2
29 //
30 // Generate first order grid function on 8 processes
31 // mpirun -np 8 ./nodal-transfer -rs 2 -rp 2 -gd 1 -o 1 -m ../../data/star.mesh
32 // Read the generated data on 4 processes and coarser mesh
33 // mpirun -np 4 ./nodal-transfer -rs 2 -rp 0 -gd 0 -snp 8 -o 1 -m ../../data/star.mesh
34 
35 #include <mfem.hpp>
36 #include <fstream>
37 #include <iostream>
38 #include <cmath>
39 #include "../common/mfem-common.hpp"
40 
41 using namespace mfem;
42 
43 class TestCoeff : public Coefficient
44 {
45 public:
46  TestCoeff() {}
47 
48  virtual
49  double Eval(ElementTransformation &T,
50  const IntegrationPoint &ip)
51  {
52  if (T.GetSpaceDim()==3)
53  {
54  double x[3];
55  Vector transip(x, 3);
56  T.Transform(ip, transip);
57  return std::sin(x[0])*std::cos(x[1]) +
58  std::sin(x[1])*std::cos(x[2]) +
59  std::sin(x[2])*std::cos(x[0]);
60  }
61  else if (T.GetSpaceDim()==2)
62  {
63  double x[2];
64  Vector transip(x, 2);
65  T.Transform(ip, transip);
66  return std::sin(x[0])*std::cos(x[1]) +
67  std::sin(x[1])*std::cos(x[0]);
68  }
69  else
70  {
71  double x;
72  Vector transip(&x,1);
73  T.Transform(ip, transip);
74  return std::sin(x)+std::cos(x);
75  }
76  }
77 };
78 
79 int main(int argc, char* argv[])
80 {
81  // Initialize MPI.
82  Mpi::Init(argc, argv);
83  int myrank = Mpi::WorldRank();
84 
85  // Parse command-line options
86  const char *mesh_file = "../../data/beam-tet.mesh";
87  int ser_ref_levels = 3;
88  int par_ref_levels = 1;
89  int order = 1;
90  int gen_data = 1;
91  int src_num_procs = 4;
92  bool visualization = true;
93 
94  OptionsParser args(argc, argv);
95  args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
96  args.AddOption(&ser_ref_levels,
97  "-rs",
98  "--refine-serial",
99  "Number of times to refine the mesh uniformly in serial.");
100  args.AddOption(&par_ref_levels,
101  "-rp",
102  "--refine-parallel",
103  "Number of times to refine the mesh uniformly in parallel.");
104  args.AddOption(&gen_data,
105  "-gd",
106  "--generate-data",
107  "Generate input data for the transfer.");
108  args.AddOption(&order,
109  "-o",
110  "--order",
111  "Order (degree) of the finite elements.");
112  args.AddOption(&src_num_procs,
113  "-snp",
114  "--src_num_procs",
115  "Number of processes for the src grid function.");
116  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
117  "--no-visualization",
118  "Enable or disable ParaView visualization.");
119  args.Parse();
120  if (!args.Good())
121  {
122  if (myrank == 0)
123  {
124  args.PrintUsage(std::cout);
125  }
126  return 1;
127  }
128 
129  if (myrank == 0)
130  {
131  args.PrintOptions(std::cout);
132  }
133 
134  // Read the (serial) mesh from the given mesh file on all processors. We
135  // can handle triangular, quadrilateral, tetrahedral and hexahedral meshes
136  // with the same code.
137  Mesh mesh(mesh_file, 1, 1);
138  int dim = mesh.SpaceDimension();
139 
140  // Refine the mesh in serial to increase the resolution. In this example
141  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
142  // a command-line parameter.
143  for (int lev = 0; lev < ser_ref_levels; lev++)
144  {
145  mesh.UniformRefinement();
146  }
147 
148  // Define a parallel mesh by a partitioning of the serial mesh. Refine
149  // this mesh further in parallel to increase the resolution. Once the
150  // parallel mesh is defined, the serial mesh can be deleted.
151  ParMesh pmesh(MPI_COMM_WORLD, mesh);
152  for (int lev = 0; lev < par_ref_levels; lev++)
153  {
154  pmesh.UniformRefinement();
155  }
156 
157  // Define the finite element spaces for the solution
158  H1_FECollection fec(order, dim);
159  ParFiniteElementSpace fespace(&pmesh, &fec, 1, Ordering::byVDIM);
160  HYPRE_Int glob_size = fespace.GlobalTrueVSize();
161  if (myrank == 0)
162  {
163  std::cout << "Number of finite element unknowns: " << glob_size
164  << std::endl;
165  }
166 
167  ParGridFunction x(&fespace); x=0.0;
168  TestCoeff prco;
169  if (gen_data)
170  {
171  Coefficient* coef[2]; coef[0]=&prco; coef[1]=&prco;
172  x.ProjectCoefficient(coef);
173 
174  // Save the grid function
175  {
176  // Save the mesh and the data
177  std::ostringstream oss;
178  oss << std::setw(10) << std::setfill('0') << myrank;
179  std::string mname="mesh_"+oss.str()+".msh";
180  std::string gname="gridfunc_"+oss.str()+".gf";
181  std::ofstream sout;
182 
183  // Save the mesh
184  sout.open(mname.c_str(),std::ios::out);
185  sout.precision(20);
186  pmesh.ParPrint(sout);
187  sout.close();
188 
189  // Save the grid function data
190  sout.open(gname.c_str(),std::ios::out);
191  sout.precision(20);
192  x.Save(sout);
193  sout.close();
194  }
195  }
196  else
197  {
198  // Read the grid function written to files and map it to the current
199  // partition scheme.
200  // x grid function will be the target of the transfer
201  // y will be utilized later for comparison
202  ParGridFunction y(&fespace);
203  Coefficient* coef[2]; coef[0]=&prco; coef[1]=&prco;
204  y.ProjectCoefficient(coef);
205 
206  // Map the src grid function
207  {
208  std::ifstream in;
210  if (dim==2)
211  {
212  map = new KDTreeNodalProjection<2>(x);
213  }
214  else
215  {
216  map = new KDTreeNodalProjection<3>(x);
217  }
218  for (int p=0; p<src_num_procs; p++)
219  {
220  std::ostringstream oss;
221  oss << std::setw(10) << std::setfill('0') << p;
222  std::string mname="mesh_"+oss.str()+".msh";
223  std::string gname="gridfunc_"+oss.str()+".gf";
224 
225  // Read the mesh
226  Mesh lmesh;
227  in.open(mname.c_str(),std::ios::in);
228  lmesh.Load(in);
229  in.close();
230 
231  in.open(gname.c_str(),std::ios::in);
232  GridFunction gf(&lmesh,in);
233  in.close();
234 
235  // Project the grid function
236  map->Project(gf,1e-8);
237  }
238  delete map;
239  }
240 
241  // Write the result into a ParaView file
242  if (visualization)
243  {
244  ParaViewDataCollection paraview_dc("GridFunc", &pmesh);
245  paraview_dc.SetPrefixPath("ParaView");
246  paraview_dc.SetLevelsOfDetail(order);
248  paraview_dc.SetCycle(0);
249  paraview_dc.SetTime(0.0);
250  paraview_dc.RegisterField("x",&x);
251  paraview_dc.RegisterField("y",&y);
252  paraview_dc.Save();
253  }
254 
255  // Compare the results
256  Vector tmpv = x;
257  tmpv -= y;
258  double l2err = mfem::InnerProduct(MPI_COMM_WORLD,tmpv,tmpv);
259  if (myrank==0)
260  {
261  std::cout<<"|l2 error|="<<sqrt(l2err)<<std::endl;
262  }
263  }
264 
265  return 0;
266 }
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetDataFormat(VTKFormat fmt)
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
Helper class for ParaView visualization data.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
void ParPrint(std::ostream &out) const
Save the mesh in a parallel mesh format.
Definition: pmesh.cpp:6295
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
virtual void Project(const Vector &coords, const Vector &src, int ordering, double lerr)=0
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
int main(int argc, char *argv[])
virtual void Load(std::istream &input, int generate_edges=0, int refine=1, bool fix_orientation=true)
Definition: mesh.hpp:666
static void Init()
Singleton creation with Mpi::Init();.
void SetTime(double t)
Set physical time (for time-dependent simulations)
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Base class for KDTreeNodalProjection.
Definition: kdtree.hpp:22
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition: mesh.hpp:1023
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:909
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:414
Class for integration point with weight.
Definition: intrules.hpp:31
int dim
Definition: ex24.cpp:53
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
void SetLevelsOfDetail(int levels_of_detail_)
Vector data type.
Definition: vector.hpp:58
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
virtual void Save() override
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.