MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
gslib.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 #include "gslib.hpp"
13 
14 #ifdef MFEM_USE_GSLIB
15 
16 // Ignore warnings from the gslib header (GCC version)
17 #ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
18 #pragma GCC diagnostic push
19 #pragma GCC diagnostic ignored "-Wunused-function"
20 #endif
21 
22 // External GSLIB header (the MFEM header is gslib.hpp)
23 namespace gslib
24 {
25 #include "gslib.h"
26 }
27 
28 #ifdef MFEM_HAVE_GCC_PRAGMA_DIAGNOSTIC
29 #pragma GCC diagnostic pop
30 #endif
31 
32 namespace mfem
33 {
34 
36  : mesh(NULL), meshsplit(NULL), ir_simplex(NULL),
37  fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
38  dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
39  avgtype(AvgType::ARITHMETIC)
40 {
41  gsl_comm = new gslib::comm;
42  cr = new gslib::crystal;
43 #ifdef MFEM_USE_MPI
44  int initialized;
45  MPI_Initialized(&initialized);
46  if (!initialized) { MPI_Init(NULL, NULL); }
47  MPI_Comm comm = MPI_COMM_WORLD;
48  comm_init(gsl_comm, comm);
49 #else
50  comm_init(gsl_comm, 0);
51 #endif
52 }
53 
55 {
56  delete gsl_comm;
57  delete cr;
58  delete ir_simplex;
59  delete meshsplit;
60 }
61 
62 #ifdef MFEM_USE_MPI
64  : mesh(NULL), meshsplit(NULL), ir_simplex(NULL),
65  fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
66  dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
67  avgtype(AvgType::ARITHMETIC)
68 {
69  gsl_comm = new gslib::comm;
70  cr = new gslib::crystal;
71  comm_init(gsl_comm, comm_);
72 }
73 #endif
74 
75 void FindPointsGSLIB::Setup(Mesh &m, const double bb_t, const double newt_tol,
76  const int npt_max)
77 {
78  MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
79  MFEM_VERIFY(m.GetNumGeometries(m.Dimension()) == 1,
80  "Mixed meshes are not currently supported in FindPointsGSLIB.");
81  MFEM_VERIFY(!(m.GetNodes()->FESpace()->IsVariableOrder()),
82  "Variable order mesh is not currently supported.");
83 
84  // call FreeData if FindPointsGSLIB::Setup has been called already
85  if (setupflag) { FreeData(); }
86 
87  crystal_init(cr, gsl_comm);
88  mesh = &m;
89  dim = mesh->Dimension();
90  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
91  unsigned dof1D = fe->GetOrder() + 1;
92  const int gt = fe->GetGeomType();
93 
94  if (gt == Geometry::TRIANGLE || gt == Geometry::TETRAHEDRON ||
95  gt == Geometry::PRISM)
96  {
98  }
99  else if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
100  {
102  }
103  else
104  {
105  MFEM_ABORT("Element type not currently supported in FindPointsGSLIB.");
106  }
107 
108  const int pts_cnt = gsl_mesh.Size()/dim,
109  NEtot = pts_cnt/(int)pow(dof1D, dim);
110 
111  if (dim == 2)
112  {
113  unsigned nr[2] = { dof1D, dof1D };
114  unsigned mr[2] = { 2*dof1D, 2*dof1D };
115  double * const elx[2] = { &gsl_mesh(0), &gsl_mesh(pts_cnt) };
116  fdata2D = findpts_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
117  pts_cnt, pts_cnt, npt_max, newt_tol);
118  }
119  else
120  {
121  unsigned nr[3] = { dof1D, dof1D, dof1D };
122  unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
123  double * const elx[3] =
124  { &gsl_mesh(0), &gsl_mesh(pts_cnt), &gsl_mesh(2*pts_cnt) };
125  fdata3D = findpts_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
126  pts_cnt, pts_cnt, npt_max, newt_tol);
127  }
128  setupflag = true;
129 }
130 
131 void FindPointsGSLIB::FindPoints(const Vector &point_pos)
132 {
133  MFEM_VERIFY(setupflag, "Use FindPointsGSLIB::Setup before finding points.");
134  points_cnt = point_pos.Size() / dim;
140 
141  if (dim == 2)
142  {
143  const double *xv_base[2];
144  xv_base[0] = point_pos.GetData();
145  xv_base[1] = point_pos.GetData() + points_cnt;
146  unsigned xv_stride[2];
147  xv_stride[0] = sizeof(double);
148  xv_stride[1] = sizeof(double);
149  findpts_2(gsl_code.GetData(), sizeof(unsigned int),
150  gsl_proc.GetData(), sizeof(unsigned int),
151  gsl_elem.GetData(), sizeof(unsigned int),
152  gsl_ref.GetData(), sizeof(double) * dim,
153  gsl_dist.GetData(), sizeof(double),
154  xv_base, xv_stride, points_cnt, fdata2D);
155  }
156  else
157  {
158  const double *xv_base[3];
159  xv_base[0] = point_pos.GetData();
160  xv_base[1] = point_pos.GetData() + points_cnt;
161  xv_base[2] = point_pos.GetData() + 2*points_cnt;
162  unsigned xv_stride[3];
163  xv_stride[0] = sizeof(double);
164  xv_stride[1] = sizeof(double);
165  xv_stride[2] = sizeof(double);
166  findpts_3(gsl_code.GetData(), sizeof(unsigned int),
167  gsl_proc.GetData(), sizeof(unsigned int),
168  gsl_elem.GetData(), sizeof(unsigned int),
169  gsl_ref.GetData(), sizeof(double) * dim,
170  gsl_dist.GetData(), sizeof(double),
171  xv_base, xv_stride, points_cnt, fdata3D);
172  }
173 
174  // Set the element number and reference position to 0 for points not found
175  for (int i = 0; i < points_cnt; i++)
176  {
177  if (gsl_code[i] == 2)
178  {
179  gsl_elem[i] = 0;
180  for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
181  }
182  }
183 
184  // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for
185  // both simplices and quads.
187 }
188 
189 void FindPointsGSLIB::FindPoints(Mesh &m, const Vector &point_pos,
190  const double bb_t, const double newt_tol,
191  const int npt_max)
192 {
193  if (!setupflag || (mesh != &m) )
194  {
195  Setup(m, bb_t, newt_tol, npt_max);
196  }
197  FindPoints(point_pos);
198 }
199 
200 void FindPointsGSLIB::Interpolate(const Vector &point_pos,
201  const GridFunction &field_in, Vector &field_out)
202 {
203  FindPoints(point_pos);
204  Interpolate(field_in, field_out);
205 }
206 
207 void FindPointsGSLIB::Interpolate(Mesh &m, const Vector &point_pos,
208  const GridFunction &field_in, Vector &field_out)
209 {
210  FindPoints(m, point_pos);
211  Interpolate(field_in, field_out);
212 }
213 
215 {
216  if (!setupflag) { return; }
217  crystal_free(cr);
218  if (dim == 2)
219  {
220  findpts_free_2(fdata2D);
221  }
222  else
223  {
224  findpts_free_3(fdata3D);
225  }
229  gsl_mesh.Destroy();
230  gsl_ref.Destroy();
231  gsl_dist.Destroy();
232  setupflag = false;
233 }
234 
236  Vector &node_vals)
237 {
238  MFEM_ASSERT(gf_in.FESpace()->GetVDim() == 1, "Scalar function expected.");
239 
240  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
241  const Geometry::Type gt = fe->GetGeomType();
242  const int NE = mesh->GetNE();
243 
244  if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
245  {
246  const GridFunction *nodes = mesh->GetNodes();
247  const FiniteElementSpace *fes = nodes->FESpace();
248  const IntegrationRule &ir = fes->GetFE(0)->GetNodes();
249  const int dof_cnt = ir.GetNPoints();
250 
251  node_vals.SetSize(NE * dof_cnt);
252 
253  const TensorBasisElement *tbe =
254  dynamic_cast<const TensorBasisElement *>(fes->GetFE(0));
255  MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
256  const Array<int> &dof_map = tbe->GetDofMap();
257 
258  int pt_id = 0;
259  Vector vals_el;
260  for (int i = 0; i < NE; i++)
261  {
262  gf_in.GetValues(i, ir, vals_el);
263  for (int j = 0; j < dof_cnt; j++)
264  {
265  node_vals(pt_id++) = vals_el(dof_map[j]);
266  }
267  }
268  }
269  else if (gt == Geometry::TRIANGLE || gt == Geometry::TETRAHEDRON ||
270  gt == Geometry::PRISM)
271  {
272  const int dof_cnt = ir_simplex->GetNPoints();
273  node_vals.SetSize(NE * dof_cnt);
274 
275  int pt_id = 0;
276  Vector vals_el;
277  for (int j = 0; j < NE; j++)
278  {
279  gf_in.GetValues(j, *ir_simplex, vals_el);
280  for (int i = 0; i < dof_cnt; i++)
281  {
282  node_vals(pt_id++) = vals_el(i);
283  }
284  }
285  }
286  else
287  {
288  MFEM_ABORT("Element type not currently supported.");
289  }
290 }
291 
293 {
294  const GridFunction *nodes = mesh->GetNodes();
295  const FiniteElementSpace *fes = nodes->FESpace();
296 
297  const int NE = mesh->GetNE(),
298  dof_cnt = fes->GetFE(0)->GetDof(),
299  pts_cnt = NE * dof_cnt;
300  gsl_mesh.SetSize(dim * pts_cnt);
301 
302  const TensorBasisElement *tbe =
303  dynamic_cast<const TensorBasisElement *>(fes->GetFE(0));
304  MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
305  Array<int> dof_map(dof_cnt);
306  const Array<int> &dm = tbe->GetDofMap();
307  if (dm.Size() > 0) { dof_map = dm; }
308  else { for (int i = 0; i < dof_cnt; i++) { dof_map[i] = i; } }
309 
310  DenseMatrix pos(dof_cnt, dim);
311  Vector posV(pos.Data(), dof_cnt * dim);
312  Array<int> xdofs(dof_cnt * dim);
313 
314  int pt_id = 0;
315  for (int i = 0; i < NE; i++)
316  {
317  fes->GetElementVDofs(i, xdofs);
318  nodes->GetSubVector(xdofs, posV);
319  for (int j = 0; j < dof_cnt; j++)
320  {
321  for (int d = 0; d < dim; d++)
322  {
323  gsl_mesh(pts_cnt * d + pt_id) = pos(dof_map[j], d);
324  }
325  pt_id++;
326  }
327  }
328 }
329 
331 {
332  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
333  const Geometry::Type gt = fe->GetGeomType();
334  const GridFunction *nodes = mesh->GetNodes();
335  const int NE = mesh->GetNE();
336  int NEsplit = 0;
337 
338  // Split the reference element into a reference submesh of quads or hexes.
339  if (gt == Geometry::TRIANGLE)
340  {
341  int Nvert = 7;
342  NEsplit = 3;
343  meshsplit = new Mesh(2, Nvert, NEsplit, 0, 2);
344 
345  const double quad_v[7][2] =
346  {
347  {0, 0}, {0.5, 0}, {1, 0}, {0, 0.5},
348  {1./3., 1./3.}, {0.5, 0.5}, {0, 1}
349  };
350  const int quad_e[3][4] =
351  {
352  {3, 4, 1, 0}, {4, 5, 2, 1}, {6, 5, 4, 3}
353  };
354 
355  for (int j = 0; j < Nvert; j++)
356  {
357  meshsplit->AddVertex(quad_v[j]);
358  }
359  for (int j = 0; j < NEsplit; j++)
360  {
361  int attribute = j + 1;
362  meshsplit->AddQuad(quad_e[j], attribute);
363  }
364  meshsplit->FinalizeQuadMesh(1, 1, true);
365  }
366  else if (gt == Geometry::TETRAHEDRON)
367  {
368  int Nvert = 15;
369  NEsplit = 4;
370  meshsplit = new Mesh(3, Nvert, NEsplit, 0, 3);
371 
372  const double hex_v[15][3] =
373  {
374  {0, 0, 0.}, {1, 0., 0.}, {0., 1., 0.}, {0, 0., 1.},
375  {0.5, 0., 0.}, {0.5, 0.5, 0.}, {0., 0.5, 0.},
376  {0., 0., 0.5}, {0.5, 0., 0.5}, {0., 0.5, 0.5},
377  {1./3., 0., 1./3.}, {1./3., 1./3., 1./3.}, {0, 1./3., 1./3.},
378  {1./3., 1./3., 0}, {0.25, 0.25, 0.25}
379  };
380  const int hex_e[4][8] =
381  {
382  {0, 4, 10, 7, 6, 13, 14, 12},
383  {4, 1, 8, 10, 13, 5, 11, 14},
384  {13, 5, 11, 14, 6, 2, 9, 12},
385  {10, 8, 3, 7, 14, 11, 9, 12}
386  };
387 
388  for (int j = 0; j < Nvert; j++)
389  {
390  meshsplit->AddVertex(hex_v[j]);
391  }
392  for (int j = 0; j < NEsplit; j++)
393  {
394  int attribute = j + 1;
395  meshsplit->AddHex(hex_e[j], attribute);
396  }
397  meshsplit->FinalizeHexMesh(1, 1, true);
398  }
399  else if (gt == Geometry::PRISM)
400  {
401  int Nvert = 14;
402  NEsplit = 3;
403  meshsplit = new Mesh(3, Nvert, NEsplit, 0, 3);
404 
405  const double hex_v[14][3] =
406  {
407  {0, 0, 0}, {0.5, 0, 0}, {1, 0, 0}, {0, 0.5, 0},
408  {1./3., 1./3., 0}, {0.5, 0.5, 0}, {0, 1, 0},
409  {0, 0, 1}, {0.5, 0, 1}, {1, 0, 1}, {0, 0.5, 1},
410  {1./3., 1./3., 1}, {0.5, 0.5, 1}, {0, 1, 1}
411  };
412  const int hex_e[3][8] =
413  {
414  {3, 4, 1, 0, 10, 11, 8, 7},
415  {4, 5, 2, 1, 11, 12, 9, 8},
416  {6, 5, 4, 3, 13, 12, 11, 10}
417  };
418 
419  for (int j = 0; j < Nvert; j++)
420  {
421  meshsplit->AddVertex(hex_v[j]);
422  }
423  for (int j = 0; j < NEsplit; j++)
424  {
425  int attribute = j + 1;
426  meshsplit->AddHex(hex_e[j], attribute);
427  }
428  meshsplit->FinalizeHexMesh(1, 1, true);
429  }
430  else { MFEM_ABORT("Unsupported geometry type."); }
431 
432  // Curve the reference submesh.
433  H1_FECollection fec(fe->GetOrder(), dim);
434  FiniteElementSpace nodal_fes(meshsplit, &fec, dim);
435  meshsplit->SetNodalFESpace(&nodal_fes);
436 
437  const int dof_cnt = nodal_fes.GetFE(0)->GetDof(),
438  pts_cnt = NEsplit * dof_cnt;
439  Vector irlist(dim * pts_cnt);
440 
441  const TensorBasisElement *tbe =
442  dynamic_cast<const TensorBasisElement *>(nodal_fes.GetFE(0));
443  MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
444  const Array<int> &dof_map = tbe->GetDofMap();
445 
446  DenseMatrix pos(dof_cnt, dim);
447  Vector posV(pos.Data(), dof_cnt * dim);
448  Array<int> xdofs(dof_cnt * dim);
449 
450  // Create an IntegrationRule on the nodes of the reference submesh.
451  ir_simplex = new IntegrationRule(pts_cnt);
452  GridFunction *nodesplit = meshsplit->GetNodes();
453  int pt_id = 0;
454  for (int i = 0; i < NEsplit; i++)
455  {
456  nodal_fes.GetElementVDofs(i, xdofs);
457  nodesplit->GetSubVector(xdofs, posV);
458  for (int j = 0; j < dof_cnt; j++)
459  {
460  for (int d = 0; d < dim; d++)
461  {
462  irlist(pts_cnt * d + pt_id) = pos(dof_map[j], d);
463  }
464  ir_simplex->IntPoint(pt_id).x = irlist(pt_id);
465  ir_simplex->IntPoint(pt_id).y = irlist(pts_cnt + pt_id);
466  if (dim == 3)
467  {
468  ir_simplex->IntPoint(pt_id).z = irlist(2*pts_cnt + pt_id);
469  }
470  pt_id++;
471  }
472  }
473 
474  // Initialize gsl_mesh with the positions of the split physical elements.
475  pt_id = 0;
476  Vector locval(dim);
477  const int tot_pts_cnt = pts_cnt*NE;
478  gsl_mesh.SetSize(tot_pts_cnt*dim);
479  for (int j = 0; j < NE; j++)
480  {
481  for (int i = 0; i < dof_cnt*NEsplit; i++)
482  {
483  const IntegrationPoint &ip = ir_simplex->IntPoint(i);
484  nodes->GetVectorValue(j, ip, locval);
485  for (int d = 0; d < dim; d++)
486  {
487  gsl_mesh(tot_pts_cnt*d + pt_id) = locval(d);
488  }
489  pt_id++;
490  }
491  }
492 }
493 
495 {
498  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
499  const Geometry::Type gt = fe->GetGeomType();
500  int NEsplit = 0;
501 
502  gsl_mfem_ref -= -1.; // map [-1, 1] to
503  gsl_mfem_ref *= 0.5; // [0, 1]
504  if (gt == Geometry::SQUARE || gt == Geometry::CUBE) { return; }
505 
506  H1_FECollection feclin(1, dim);
507  FiniteElementSpace nodal_fes_lin(meshsplit, &feclin, dim);
508  GridFunction gf_lin(&nodal_fes_lin);
509 
510  if (gt == Geometry::TRIANGLE)
511  {
512  const double quad_v[7][2] =
513  {
514  {0, 0}, {0.5, 0}, {1, 0}, {0, 0.5},
515  {1./3., 1./3.}, {0.5, 0.5}, {0, 1}
516  };
517  for (int k = 0; k < dim; k++)
518  {
519  for (int j = 0; j < gf_lin.Size()/dim; j++)
520  {
521  gf_lin(j+k*gf_lin.Size()/dim) = quad_v[j][k];
522  }
523  }
524  NEsplit = 3;
525  }
526  else if (gt == Geometry::TETRAHEDRON)
527  {
528  const double hex_v[15][3] =
529  {
530  {0, 0, 0.}, {1, 0., 0.}, {0., 1., 0.}, {0, 0., 1.},
531  {0.5, 0., 0.}, {0.5, 0.5, 0.}, {0., 0.5, 0.},
532  {0., 0., 0.5}, {0.5, 0., 0.5}, {0., 0.5, 0.5},
533  {1./3., 0., 1./3.}, {1./3., 1./3., 1./3.}, {0, 1./3., 1./3.},
534  {1./3., 1./3., 0}, {0.25, 0.25, 0.25}
535  };
536  for (int k = 0; k < dim; k++)
537  {
538  for (int j = 0; j < gf_lin.Size()/dim; j++)
539  {
540  gf_lin(j+k*gf_lin.Size()/dim) = hex_v[j][k];
541  }
542  }
543  NEsplit = 4;
544  }
545  else if (gt == Geometry::PRISM)
546  {
547  const double hex_v[14][3] =
548  {
549  {0, 0, 0}, {0.5, 0, 0}, {1, 0, 0}, {0, 0.5, 0},
550  {1./3., 1./3., 0}, {0.5, 0.5, 0}, {0, 1, 0},
551  {0, 0, 1}, {0.5, 0, 1}, {1, 0, 1}, {0, 0.5, 1},
552  {1./3., 1./3., 1}, {0.5, 0.5, 1}, {0, 1, 1}
553  };
554  for (int k = 0; k < dim; k++)
555  {
556  for (int j = 0; j < gf_lin.Size()/dim; j++)
557  {
558  gf_lin(j+k*gf_lin.Size()/dim) = hex_v[j][k];
559  }
560  }
561  NEsplit = 3;
562  }
563  else
564  {
565  MFEM_ABORT("Element type not currently supported.");
566  }
567 
568  // Simplices are split into quads/hexes for GSLIB. For MFEM, we need to find
569  // the original element number and map the rst from micro to macro element.
570  for (int i = 0; i < points_cnt; i++)
571  {
572  if (gsl_code[i] == 2) { continue; }
573  int local_elem = gsl_elem[i]%NEsplit;
574  gsl_mfem_elem[i] = (gsl_elem[i] - local_elem)/NEsplit; // macro element number
575 
576  IntegrationPoint ip;
577  Vector mfem_ref(gsl_mfem_ref.GetData()+i*dim, dim);
578  ip.Set2(mfem_ref.GetData());
579  if (dim == 3) { ip.z = mfem_ref(2); }
580  gf_lin.GetVectorValue(local_elem, ip, mfem_ref); // map to rst of macro element
581  }
582 }
583 
585  Vector &field_out)
586 {
587  const int gf_order = field_in.FESpace()->GetFE(0)->GetOrder(),
588  mesh_order = mesh->GetNodalFESpace()->GetFE(0)->GetOrder();
589 
590  const FiniteElementCollection *fec_in = field_in.FESpace()->FEColl();
591  const H1_FECollection *fec_h1 = dynamic_cast<const H1_FECollection *>(fec_in);
592  const L2_FECollection *fec_l2 = dynamic_cast<const L2_FECollection *>(fec_in);
593 
594  if (fec_h1 && gf_order == mesh_order &&
595  fec_h1->GetBasisType() == BasisType::GaussLobatto &&
596  !field_in.FESpace()->IsVariableOrder())
597  {
598  InterpolateH1(field_in, field_out);
599  return;
600  }
601  else
602  {
603  InterpolateGeneral(field_in, field_out);
604  if (!fec_l2 || avgtype == AvgType::NONE) { return; }
605  }
606 
607  // For points on element borders, project the L2 GridFunction to H1 and
608  // re-interpolate.
609  if (fec_l2)
610  {
611  Array<int> indl2;
612  for (int i = 0; i < points_cnt; i++)
613  {
614  if (gsl_code[i] == 1) { indl2.Append(i); }
615  }
616  int borderPts = indl2.Size();
617 #ifdef MFEM_USE_MPI
618  MPI_Allreduce(MPI_IN_PLACE, &borderPts, 1, MPI_INT, MPI_SUM, gsl_comm->c);
619 #endif
620  if (borderPts == 0) { return; } // no points on element borders
621 
622 
623  Vector field_out_l2(field_out.Size());
624  VectorGridFunctionCoefficient field_in_dg(&field_in);
625  int gf_order_h1 = std::max(gf_order, 1); // H1 should be at least order 1
626  H1_FECollection fec(gf_order_h1, dim);
627  const int ncomp = field_in.FESpace()->GetVDim();
628  FiniteElementSpace fes(mesh, &fec, ncomp);
629  GridFunction field_in_h1(&fes);
630 
631  if (avgtype == AvgType::ARITHMETIC)
632  {
633  field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::ARITHMETIC);
634  }
635  else if (avgtype == AvgType::HARMONIC)
636  {
637  field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::HARMONIC);
638  }
639  else
640  {
641  MFEM_ABORT("Invalid averaging type.");
642  }
643 
644  if (gf_order_h1 == mesh_order) // basis is GaussLobatto by default
645  {
646  InterpolateH1(field_in_h1, field_out_l2);
647  }
648  else
649  {
650  InterpolateGeneral(field_in_h1, field_out_l2);
651  }
652 
653  // Copy interpolated values for the points on element border
654  for (int j = 0; j < ncomp; j++)
655  {
656  for (int i = 0; i < indl2.Size(); i++)
657  {
658  int idx = indl2[i] + j*points_cnt;
659  field_out(idx) = field_out_l2(idx);
660  }
661  }
662  }
663 }
664 
666  Vector &field_out)
667 {
668  FiniteElementSpace ind_fes(mesh, field_in.FESpace()->FEColl());
669  GridFunction field_in_scalar(&ind_fes);
670  Vector node_vals;
671 
672  const int ncomp = field_in.FESpace()->GetVDim(),
673  points_fld = field_in.Size() / ncomp,
675 
676  field_out.SetSize(points_cnt*ncomp);
677  field_out = default_interp_value;
678 
679  for (int i = 0; i < ncomp; i++)
680  {
681  const int dataptrin = i*points_fld,
682  dataptrout = i*points_cnt;
683  field_in_scalar.NewDataAndSize(field_in.GetData()+dataptrin, points_fld);
684  GetNodeValues(field_in_scalar, node_vals);
685 
686  if (dim==2)
687  {
688  findpts_eval_2(field_out.GetData()+dataptrout, sizeof(double),
689  gsl_code.GetData(), sizeof(unsigned int),
690  gsl_proc.GetData(), sizeof(unsigned int),
691  gsl_elem.GetData(), sizeof(unsigned int),
692  gsl_ref.GetData(), sizeof(double) * dim,
693  points_cnt, node_vals.GetData(), fdata2D);
694  }
695  else
696  {
697  findpts_eval_3(field_out.GetData()+dataptrout, sizeof(double),
698  gsl_code.GetData(), sizeof(unsigned int),
699  gsl_proc.GetData(), sizeof(unsigned int),
700  gsl_elem.GetData(), sizeof(unsigned int),
701  gsl_ref.GetData(), sizeof(double) * dim,
702  points_cnt, node_vals.GetData(), fdata3D);
703  }
704  }
705 }
706 
708  Vector &field_out)
709 {
710  int ncomp = field_in.VectorDim(),
711  nptorig = points_cnt,
712  npt = points_cnt;
713 
714  field_out.SetSize(points_cnt*ncomp);
715  field_out = default_interp_value;
716 
717  if (gsl_comm->np == 1) // serial
718  {
719  for (int index = 0; index < npt; index++)
720  {
721  if (gsl_code[index] == 2) { continue; }
722  IntegrationPoint ip;
724  if (dim == 3) { ip.z = gsl_mfem_ref(index*dim + 2); }
725  Vector localval(ncomp);
726  field_in.GetVectorValue(gsl_mfem_elem[index], ip, localval);
727  for (int i = 0; i < ncomp; i++)
728  {
729  field_out(index + i*npt) = localval(i);
730  }
731  }
732  }
733  else // parallel
734  {
735  // Determine number of points to be sent
736  int nptsend = 0;
737  for (int index = 0; index < npt; index++)
738  {
739  if (gsl_code[index] != 2) { nptsend +=1; }
740  }
741 
742  // Pack data to send via crystal router
743  struct gslib::array *outpt = new gslib::array;
744  struct out_pt { double r[3], ival; uint index, el, proc; };
745  struct out_pt *pt;
746  array_init(struct out_pt, outpt, nptsend);
747  outpt->n=nptsend;
748  pt = (struct out_pt *)outpt->ptr;
749  for (int index = 0; index < npt; index++)
750  {
751  if (gsl_code[index] == 2) { continue; }
752  for (int d = 0; d < dim; ++d) { pt->r[d]= gsl_mfem_ref(index*dim + d); }
753  pt->index = index;
754  pt->proc = gsl_proc[index];
755  pt->el = gsl_mfem_elem[index];
756  ++pt;
757  }
758 
759  // Transfer data to target MPI ranks
760  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
761 
762  if (ncomp == 1)
763  {
764  // Interpolate the grid function
765  npt = outpt->n;
766  pt = (struct out_pt *)outpt->ptr;
767  for (int index = 0; index < npt; index++)
768  {
769  IntegrationPoint ip;
770  ip.Set3(&pt->r[0]);
771  pt->ival = field_in.GetValue(pt->el, ip, 1);
772  ++pt;
773  }
774 
775  // Transfer data back to source MPI rank
776  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
777  npt = outpt->n;
778  pt = (struct out_pt *)outpt->ptr;
779  for (int index = 0; index < npt; index++)
780  {
781  field_out(pt->index) = pt->ival;
782  ++pt;
783  }
784  array_free(outpt);
785  delete outpt;
786  }
787  else // ncomp > 1
788  {
789  // Interpolate data and store in a Vector
790  npt = outpt->n;
791  pt = (struct out_pt *)outpt->ptr;
792  Vector vec_int_vals(npt*ncomp);
793  for (int index = 0; index < npt; index++)
794  {
795  IntegrationPoint ip;
796  ip.Set3(&pt->r[0]);
797  Vector localval(vec_int_vals.GetData()+index*ncomp, ncomp);
798  field_in.GetVectorValue(pt->el, ip, localval);
799  ++pt;
800  }
801 
802  // Save index and proc data in a struct
803  struct gslib::array *savpt = new gslib::array;
804  struct sav_pt { uint index, proc; };
805  struct sav_pt *spt;
806  array_init(struct sav_pt, savpt, npt);
807  savpt->n=npt;
808  spt = (struct sav_pt *)savpt->ptr;
809  pt = (struct out_pt *)outpt->ptr;
810  for (int index = 0; index < npt; index++)
811  {
812  spt->index = pt->index;
813  spt->proc = pt->proc;
814  ++pt; ++spt;
815  }
816 
817  array_free(outpt);
818  delete outpt;
819 
820  // Copy data from save struct to send struct and send component wise
821  struct gslib::array *sendpt = new gslib::array;
822  struct send_pt { double ival; uint index, proc; };
823  struct send_pt *sdpt;
824  for (int j = 0; j < ncomp; j++)
825  {
826  array_init(struct send_pt, sendpt, npt);
827  sendpt->n=npt;
828  spt = (struct sav_pt *)savpt->ptr;
829  sdpt = (struct send_pt *)sendpt->ptr;
830  for (int index = 0; index < npt; index++)
831  {
832  sdpt->index = spt->index;
833  sdpt->proc = spt->proc;
834  sdpt->ival = vec_int_vals(j + index*ncomp);
835  ++sdpt; ++spt;
836  }
837 
838  sarray_transfer(struct send_pt, sendpt, proc, 1, cr);
839  sdpt = (struct send_pt *)sendpt->ptr;
840  for (int index = 0; index < nptorig; index++)
841  {
842  int idx = sdpt->index + j*nptorig;
843  field_out(idx) = sdpt->ival;
844  ++sdpt;
845  }
846  array_free(sendpt);
847  }
848  array_free(savpt);
849  delete sendpt;
850  delete savpt;
851  } // ncomp > 1
852  } // parallel
853 }
854 
855 void OversetFindPointsGSLIB::Setup(Mesh &m, const int meshid,
856  GridFunction *gfmax,
857  const double bb_t, const double newt_tol,
858  const int npt_max)
859 {
860  MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
861  MFEM_VERIFY(m.GetNumGeometries(m.Dimension()) == 1,
862  "Mixed meshes are not currently supported in FindPointsGSLIB.");
863  MFEM_VERIFY(!(m.GetNodes()->FESpace()->IsVariableOrder()),
864  "Variable order mesh is not currently supported.");
865 
866  // FreeData if OversetFindPointsGSLIB::Setup has been called already
867  if (setupflag) { FreeData(); }
868 
869  crystal_init(cr, gsl_comm);
870  mesh = &m;
871  dim = mesh->Dimension();
872  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
873  unsigned dof1D = fe->GetOrder() + 1;
874  int gt = fe->GetGeomType();
875 
876  if (gt == Geometry::TRIANGLE || gt == Geometry::TETRAHEDRON ||
877  gt == Geometry::PRISM)
878  {
880  }
881  else if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
882  {
884  }
885  else
886  {
887  MFEM_ABORT("Element type not currently supported in FindPointsGSLIB.");
888  }
889 
890  MFEM_ASSERT(meshid>=0, " The ID should be greater than or equal to 0.");
891 
892  const int pts_cnt = gsl_mesh.Size()/dim,
893  NEtot = pts_cnt/(int)pow(dof1D, dim);
894 
895  distfint.SetSize(pts_cnt);
896  if (!gfmax)
897  {
898  distfint = 0.;
899  }
900  else
901  {
902  GetNodeValues(*gfmax, distfint);
903  }
904  u_meshid = (unsigned int)meshid;
905 
906  if (dim == 2)
907  {
908  unsigned nr[2] = { dof1D, dof1D };
909  unsigned mr[2] = { 2*dof1D, 2*dof1D };
910  double * const elx[2] = { &gsl_mesh(0), &gsl_mesh(pts_cnt) };
911  fdata2D = findptsms_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
912  pts_cnt, pts_cnt, npt_max, newt_tol,
913  &u_meshid, &distfint(0));
914  }
915  else
916  {
917  unsigned nr[3] = { dof1D, dof1D, dof1D };
918  unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
919  double * const elx[3] =
920  { &gsl_mesh(0), &gsl_mesh(pts_cnt), &gsl_mesh(2*pts_cnt) };
921  fdata3D = findptsms_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
922  pts_cnt, pts_cnt, npt_max, newt_tol,
923  &u_meshid, &distfint(0));
924  }
925  setupflag = true;
926  overset = true;
927 }
928 
930  Array<unsigned int> &point_id)
931 {
932  MFEM_VERIFY(setupflag, "Use OversetFindPointsGSLIB::Setup before "
933  "finding points.");
934  MFEM_VERIFY(overset, " Please setup FindPoints for overlapping grids.");
935  points_cnt = point_pos.Size() / dim;
936  unsigned int match = 0; // Don't find points in the mesh if point_id = mesh_id
937 
943 
944  if (dim == 2)
945  {
946  const double *xv_base[2];
947  xv_base[0] = point_pos.GetData();
948  xv_base[1] = point_pos.GetData() + points_cnt;
949  unsigned xv_stride[2];
950  xv_stride[0] = sizeof(double);
951  xv_stride[1] = sizeof(double);
952  findptsms_2(gsl_code.GetData(), sizeof(unsigned int),
953  gsl_proc.GetData(), sizeof(unsigned int),
954  gsl_elem.GetData(), sizeof(unsigned int),
955  gsl_ref.GetData(), sizeof(double) * dim,
956  gsl_dist.GetData(), sizeof(double),
957  xv_base, xv_stride,
958  point_id.GetData(), sizeof(unsigned int), &match,
960  }
961  else
962  {
963  const double *xv_base[3];
964  xv_base[0] = point_pos.GetData();
965  xv_base[1] = point_pos.GetData() + points_cnt;
966  xv_base[2] = point_pos.GetData() + 2*points_cnt;
967  unsigned xv_stride[3];
968  xv_stride[0] = sizeof(double);
969  xv_stride[1] = sizeof(double);
970  xv_stride[2] = sizeof(double);
971  findptsms_3(gsl_code.GetData(), sizeof(unsigned int),
972  gsl_proc.GetData(), sizeof(unsigned int),
973  gsl_elem.GetData(), sizeof(unsigned int),
974  gsl_ref.GetData(), sizeof(double) * dim,
975  gsl_dist.GetData(), sizeof(double),
976  xv_base, xv_stride,
977  point_id.GetData(), sizeof(unsigned int), &match,
979  }
980 
981  // Set the element number and reference position to 0 for points not found
982  for (int i = 0; i < points_cnt; i++)
983  {
984  if (gsl_code[i] == 2)
985  {
986  gsl_elem[i] = 0;
987  for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
988  }
989  }
990 
991  // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for both
992  // simplices and quads.
994 }
995 
997  Array<unsigned int> &point_id,
998  const GridFunction &field_in,
999  Vector &field_out)
1000 {
1001  FindPoints(point_pos, point_id);
1002  Interpolate(field_in, field_out);
1003 }
1004 
1005 
1006 } // namespace mfem
1007 
1008 #endif // MFEM_USE_GSLIB
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:584
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
Array< unsigned int > gsl_elem
Definition: gslib.hpp:60
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1662
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:455
virtual ~FindPointsGSLIB()
Definition: gslib.cpp:54
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:465
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5877
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2569
Array< unsigned int > gsl_mfem_elem
Definition: gslib.hpp:60
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:534
T * GetData()
Returns the data.
Definition: array.hpp:112
unsigned int uint
Definition: gecko.hpp:204
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
int VectorDim() const
Definition: gridfunc.cpp:321
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
Array< unsigned int > gsl_proc
Definition: gslib.hpp:60
void FinalizeHexMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a hexahedral Mesh.
Definition: mesh.cpp:2841
virtual void GetNodeValues(const GridFunction &gf_in, Vector &node_vals)
Get GridFunction from MFEM format to GSLIB format.
Definition: gslib.cpp:235
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
virtual void MapRefPosAndElemIndices()
Definition: gslib.cpp:494
struct gslib::findpts_data_3 * fdata3D
Definition: gslib.hpp:56
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1601
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
virtual void FreeData()
Definition: gslib.cpp:214
int GetBasisType() const
Definition: fe_coll.hpp:239
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:75
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
double default_interp_value
Definition: gslib.hpp:63
struct gslib::crystal * cr
Definition: gslib.hpp:57
void Set2(const double x1, const double x2)
Definition: intrules.hpp:80
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:510
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:557
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
int Dimension() const
Definition: mesh.hpp:999
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:70
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id)
Definition: gslib.cpp:929
Array< unsigned int > gsl_code
Definition: gslib.hpp:60
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
virtual void GetSimplexNodalCoordinates()
Definition: gslib.cpp:330
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
void FindPoints(const Vector &point_pos)
Definition: gslib.cpp:131
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
int AddHex(int v1, int v2, int v3, int v4, int v5, int v6, int v7, int v8, int attr=1)
Definition: mesh.cpp:1725
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:707
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
virtual void GetQuadHexNodalCoordinates()
Definition: gslib.cpp:292
struct gslib::findpts_data_2 * fdata2D
Definition: gslib.hpp:55
const Array< int > & GetDofMap() const
Get an Array&lt;int&gt; that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe_base.hpp:1181
Class for integration point with weight.
Definition: intrules.hpp:25
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5307
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1956
int dim
Definition: ex24.cpp:53
void Destroy()
Destroy a vector.
Definition: vector.hpp:598
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2783
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
struct gslib::comm * gsl_comm
Definition: gslib.hpp:58
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:996
IntegrationRule * ir_simplex
Definition: gslib.hpp:54
Vector data type.
Definition: vector.hpp:60
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:577
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
Vector coefficient defined by a vector GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void Setup(Mesh &m, const int meshid, GridFunction *gfmax=NULL, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:855
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition: gslib.cpp:665
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:438
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:284