MFEM  v4.5.2
Finite element discretization library
gslib.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 #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),
37  fec_map_lin(NULL),
38  fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
39  dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
40  avgtype(AvgType::ARITHMETIC), bdr_tol(1e-8)
41 {
42  mesh_split.SetSize(4);
43  ir_split.SetSize(4);
45  gf_rst_map.SetSize(4);
46  for (int i = 0; i < mesh_split.Size(); i++)
47  {
48  mesh_split[i] = NULL;
49  ir_split[i] = NULL;
50  fes_rst_map[i] = NULL;
51  gf_rst_map[i] = NULL;
52  }
53 
54  gsl_comm = new gslib::comm;
55  cr = new gslib::crystal;
56 #ifdef MFEM_USE_MPI
57  int initialized;
58  MPI_Initialized(&initialized);
59  if (!initialized) { MPI_Init(NULL, NULL); }
60  MPI_Comm comm = MPI_COMM_WORLD;
61  comm_init(gsl_comm, comm);
62 #else
63  comm_init(gsl_comm, 0);
64 #endif
65 }
66 
68 {
69  delete gsl_comm;
70  delete cr;
71  for (int i = 0; i < 4; i++)
72  {
73  if (mesh_split[i]) { delete mesh_split[i]; mesh_split[i] = NULL; }
74  if (ir_split[i]) { delete ir_split[i]; ir_split[i] = NULL; }
75  if (fes_rst_map[i]) { delete fes_rst_map[i]; fes_rst_map[i] = NULL; }
76  if (gf_rst_map[i]) { delete gf_rst_map[i]; gf_rst_map[i] = NULL; }
77  }
78  if (fec_map_lin) { delete fec_map_lin; fec_map_lin = NULL; }
79 }
80 
81 #ifdef MFEM_USE_MPI
83  : mesh(NULL),
84  fec_map_lin(NULL),
85  fdata2D(NULL), fdata3D(NULL), cr(NULL), gsl_comm(NULL),
86  dim(-1), points_cnt(0), setupflag(false), default_interp_value(0),
87  avgtype(AvgType::ARITHMETIC), bdr_tol(1e-8)
88 {
89  mesh_split.SetSize(4);
90  ir_split.SetSize(4);
92  gf_rst_map.SetSize(4);
93  for (int i = 0; i < mesh_split.Size(); i++)
94  {
95  mesh_split[i] = NULL;
96  ir_split[i] = NULL;
97  fes_rst_map[i] = NULL;
98  gf_rst_map[i] = NULL;
99  }
100 
101  gsl_comm = new gslib::comm;
102  cr = new gslib::crystal;
103  comm_init(gsl_comm, comm_);
104 }
105 #endif
106 
107 void FindPointsGSLIB::Setup(Mesh &m, const double bb_t, const double newt_tol,
108  const int npt_max)
109 {
110  MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
111  MFEM_VERIFY(!(m.GetNodes()->FESpace()->IsVariableOrder()),
112  "Variable order mesh is not currently supported.");
113 
114  // call FreeData if FindPointsGSLIB::Setup has been called already
115  if (setupflag) { FreeData(); }
116 
117  crystal_init(cr, gsl_comm);
118  mesh = &m;
119  dim = mesh->Dimension();
120  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
121  unsigned dof1D = fe->GetOrder() + 1;
122 
124  if (dim == 2)
125  {
126  if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
127  ir_split[0] = new IntegrationRule(3*pow(dof1D, dim));
129  }
130  else if (dim == 3)
131  {
132 
133  if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
134  ir_split[1] = new IntegrationRule(4*pow(dof1D, dim));
136 
137  if (ir_split[2]) { delete ir_split[2]; ir_split[2] = NULL; }
138  ir_split[2] = new IntegrationRule(3*pow(dof1D, dim));
140 
141  if (ir_split[3]) { delete ir_split[3]; ir_split[3] = NULL; }
142  ir_split[3] = new IntegrationRule(8*pow(dof1D, dim));
144  }
145 
147 
148  const int pts_cnt = gsl_mesh.Size()/dim,
149  NEtot = pts_cnt/(int)pow(dof1D, dim);
150 
151  if (dim == 2)
152  {
153  unsigned nr[2] = { dof1D, dof1D };
154  unsigned mr[2] = { 2*dof1D, 2*dof1D };
155  double * const elx[2] = { &gsl_mesh(0), &gsl_mesh(pts_cnt) };
156  fdata2D = findpts_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
157  pts_cnt, pts_cnt, npt_max, newt_tol);
158  }
159  else
160  {
161  unsigned nr[3] = { dof1D, dof1D, dof1D };
162  unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
163  double * const elx[3] =
164  { &gsl_mesh(0), &gsl_mesh(pts_cnt), &gsl_mesh(2*pts_cnt) };
165  fdata3D = findpts_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
166  pts_cnt, pts_cnt, npt_max, newt_tol);
167  }
168  setupflag = true;
169 }
170 
171 void FindPointsGSLIB::FindPoints(const Vector &point_pos,
172  int point_pos_ordering)
173 {
174  MFEM_VERIFY(setupflag, "Use FindPointsGSLIB::Setup before finding points.");
175  points_cnt = point_pos.Size() / dim;
181 
182  auto xvFill = [&](const double *xv_base[], unsigned xv_stride[], int dim)
183  {
184  for (int d = 0; d < dim; d++)
185  {
186  if (point_pos_ordering == Ordering::byNODES)
187  {
188  xv_base[d] = point_pos.GetData() + d*points_cnt;
189  xv_stride[d] = sizeof(double);
190  }
191  else
192  {
193  xv_base[d] = point_pos.GetData() + d;
194  xv_stride[d] = dim*sizeof(double);
195  }
196  }
197  };
198  if (dim == 2)
199  {
200  const double *xv_base[2];
201  unsigned xv_stride[2];
202  xvFill(xv_base, xv_stride, dim);
203  findpts_2(gsl_code.GetData(), sizeof(unsigned int),
204  gsl_proc.GetData(), sizeof(unsigned int),
205  gsl_elem.GetData(), sizeof(unsigned int),
206  gsl_ref.GetData(), sizeof(double) * dim,
207  gsl_dist.GetData(), sizeof(double),
208  xv_base, xv_stride, points_cnt, fdata2D);
209  }
210  else // dim == 3
211  {
212  const double *xv_base[3];
213  unsigned xv_stride[3];
214  xvFill(xv_base, xv_stride, dim);
215  findpts_3(gsl_code.GetData(), sizeof(unsigned int),
216  gsl_proc.GetData(), sizeof(unsigned int),
217  gsl_elem.GetData(), sizeof(unsigned int),
218  gsl_ref.GetData(), sizeof(double) * dim,
219  gsl_dist.GetData(), sizeof(double),
220  xv_base, xv_stride, points_cnt, fdata3D);
221  }
222 
223  // Set the element number and reference position to 0 for points not found
224  for (int i = 0; i < points_cnt; i++)
225  {
226  if (gsl_code[i] == 2 ||
227  (gsl_code[i] == 1 && gsl_dist(i) > bdr_tol))
228  {
229  gsl_elem[i] = 0;
230  for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
231  gsl_code[i] = 2;
232  }
233  }
234 
235  // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for
236  // both simplices and quads.
238 }
239 
240 void FindPointsGSLIB::FindPoints(Mesh &m, const Vector &point_pos,
241  int point_pos_ordering, const double bb_t,
242  const double newt_tol, const int npt_max)
243 {
244  if (!setupflag || (mesh != &m) )
245  {
246  Setup(m, bb_t, newt_tol, npt_max);
247  }
248  FindPoints(point_pos, point_pos_ordering);
249 }
250 
251 void FindPointsGSLIB::Interpolate(const Vector &point_pos,
252  const GridFunction &field_in, Vector &field_out,
253  int point_pos_ordering)
254 {
255  FindPoints(point_pos, point_pos_ordering);
256  Interpolate(field_in, field_out);
257 }
258 
259 void FindPointsGSLIB::Interpolate(Mesh &m, const Vector &point_pos,
260  const GridFunction &field_in, Vector &field_out,
261  int point_pos_ordering)
262 {
263  FindPoints(m, point_pos, point_pos_ordering);
264  Interpolate(field_in, field_out);
265 }
266 
268 {
269  if (!setupflag) { return; }
270  crystal_free(cr);
271  if (dim == 2)
272  {
273  findpts_free_2(fdata2D);
274  }
275  else
276  {
277  findpts_free_3(fdata3D);
278  }
282  gsl_mesh.Destroy();
283  gsl_ref.Destroy();
284  gsl_dist.Destroy();
285  for (int i = 0; i < 4; i++)
286  {
287  if (mesh_split[i]) { delete mesh_split[i]; mesh_split[i] = NULL; }
288  if (ir_split[i]) { delete ir_split[i]; ir_split[i] = NULL; }
289  if (fes_rst_map[i]) { delete fes_rst_map[i]; fes_rst_map[i] = NULL; }
290  if (gf_rst_map[i]) { delete gf_rst_map[i]; gf_rst_map[i] = NULL; }
291  }
292  if (fec_map_lin) { delete fec_map_lin; fec_map_lin = NULL; }
293  setupflag = false;
294 }
295 
297 {
298  fec_map_lin = new H1_FECollection(1, dim);
299  if (mesh->Dimension() == 2)
300  {
301  int Nvert = 7;
302  int NEsplit = 3;
303  mesh_split[0] = new Mesh(2, Nvert, NEsplit, 0, 2);
304 
305  const double quad_v[7][2] =
306  {
307  {0, 0}, {0.5, 0}, {1, 0}, {0, 0.5},
308  {1./3., 1./3.}, {0.5, 0.5}, {0, 1}
309  };
310  const int quad_e[3][4] =
311  {
312  {0, 1, 4, 3}, {1, 2, 5, 4}, {3, 4, 5, 6}
313  };
314 
315  for (int j = 0; j < Nvert; j++)
316  {
317  mesh_split[0]->AddVertex(quad_v[j]);
318  }
319  for (int j = 0; j < NEsplit; j++)
320  {
321  int attribute = j + 1;
322  mesh_split[0]->AddQuad(quad_e[j], attribute);
323  }
324  mesh_split[0]->FinalizeQuadMesh(1, 1, true);
325 
327  gf_rst_map[0] = new GridFunction(fes_rst_map[0]);
328  const int npt = gf_rst_map[0]->Size()/dim;
329  for (int k = 0; k < dim; k++)
330  {
331  for (int j = 0; j < npt; j++)
332  {
333  (*gf_rst_map[0])(j+k*npt) = quad_v[j][k];
334  }
335  }
336  }
337  else if (mesh->Dimension() == 3)
338  {
339  // Tetrahedron
340  {
341  int Nvert = 15;
342  int NEsplit = 4;
343  mesh_split[1] = new Mesh(3, Nvert, NEsplit, 0, 3);
344 
345  const double hex_v[15][3] =
346  {
347  {0, 0, 0.}, {1, 0., 0.}, {0., 1., 0.}, {0, 0., 1.},
348  {0.5, 0., 0.}, {0.5, 0.5, 0.}, {0., 0.5, 0.},
349  {0., 0., 0.5}, {0.5, 0., 0.5}, {0., 0.5, 0.5},
350  {1./3., 0., 1./3.}, {1./3., 1./3., 1./3.}, {0, 1./3., 1./3.},
351  {1./3., 1./3., 0}, {0.25, 0.25, 0.25}
352  };
353  const int hex_e[4][8] =
354  {
355  {7, 10, 4, 0, 12, 14, 13, 6},
356  {10, 8, 1, 4, 14, 11, 5, 13},
357  {14, 11, 5, 13, 12, 9, 2, 6},
358  {7, 3, 8, 10, 12, 9, 11, 14}
359  };
360 
361  for (int j = 0; j < Nvert; j++)
362  {
363  mesh_split[1]->AddVertex(hex_v[j]);
364  }
365  for (int j = 0; j < NEsplit; j++)
366  {
367  int attribute = j + 1;
368  mesh_split[1]->AddHex(hex_e[j], attribute);
369  }
370  mesh_split[1]->FinalizeHexMesh(1, 1, true);
371 
373  gf_rst_map[1] = new GridFunction(fes_rst_map[1]);
374  const int npt = gf_rst_map[1]->Size()/dim;
375  for (int k = 0; k < dim; k++)
376  {
377  for (int j = 0; j < npt; j++)
378  {
379  (*gf_rst_map[1])(j+k*npt) = hex_v[j][k];
380  }
381  }
382  }
383  // Prism
384  {
385  int Nvert = 14;
386  int NEsplit = 3;
387  mesh_split[2] = new Mesh(3, Nvert, NEsplit, 0, 3);
388 
389  const double hex_v[14][3] =
390  {
391  {0, 0, 0}, {0.5, 0, 0}, {1, 0, 0}, {0, 0.5, 0},
392  {1./3., 1./3., 0}, {0.5, 0.5, 0}, {0, 1, 0},
393  {0, 0, 1}, {0.5, 0, 1}, {1, 0, 1}, {0, 0.5, 1},
394  {1./3., 1./3., 1}, {0.5, 0.5, 1}, {0, 1, 1}
395  };
396  const int hex_e[3][8] =
397  {
398  {0, 1, 4, 3, 7, 8, 11, 10},
399  {1, 2, 5, 4, 8, 9, 12, 11},
400  {3, 4, 5, 6, 10, 11, 12, 13}
401  };
402 
403  for (int j = 0; j < Nvert; j++)
404  {
405  mesh_split[2]->AddVertex(hex_v[j]);
406  }
407  for (int j = 0; j < NEsplit; j++)
408  {
409  int attribute = j + 1;
410  mesh_split[2]->AddHex(hex_e[j], attribute);
411  }
412  mesh_split[2]->FinalizeHexMesh(1, 1, true);
413 
415  gf_rst_map[2] = new GridFunction(fes_rst_map[2]);
416  const int npt = gf_rst_map[2]->Size()/dim;
417  for (int k = 0; k < dim; k++)
418  {
419  for (int j = 0; j < npt; j++)
420  {
421  (*gf_rst_map[2])(j+k*npt) = hex_v[j][k];
422  }
423  }
424  }
425  // Pyramid
426  {
427  int Nvert = 23;
428  int NEsplit = 8;
429  mesh_split[3] = new Mesh(3, Nvert, NEsplit, 0, 3);
430 
431  const double hex_v[23][3] =
432  {
433  {0.0000, 0.0000, 0.0000}, {0.5000, 0.0000, 0.0000},
434  {0.0000, 0.0000, 0.5000}, {0.3333, 0.0000, 0.3333},
435  {0.0000, 0.5000, 0.0000}, {0.3333, 0.3333, 0.0000},
436  {0.0000, 0.3333, 0.3333}, {0.2500, 0.2500, 0.2500},
437  {1.0000, 0.0000, 0.0000}, {0.5000, 0.0000, 0.5000},
438  {0.5000, 0.5000, 0.0000}, {0.3333, 0.3333, 0.3333},
439  {0.0000, 1.0000, 0.0000}, {0.0000, 0.5000, 0.5000},
440  {0.0000, 0.0000, 1.0000}, {1.0000, 0.5000, 0.0000},
441  {0.6667, 0.3333, 0.3333}, {0.6667, 0.6667, 0.0000},
442  {0.5000, 0.5000, 0.2500}, {1.0000, 1.0000, 0.0000},
443  {0.5000, 0.5000, 0.5000}, {0.5000, 1.0000, 0.0000},
444  {0.3333, 0.6667, 0.3333}
445  };
446  const int hex_e[8][8] =
447  {
448  {2, 3, 1, 0, 6, 7, 5, 4}, {3, 9, 8, 1, 7, 11, 10, 5},
449  {7, 11, 10, 5, 6, 13, 12, 4}, {2, 14, 9, 3, 6, 13, 11, 7},
450  {9, 16, 15, 8, 11, 18, 17, 10}, {16, 20, 19, 15, 18, 22, 21, 17},
451  {18, 22, 21, 17, 11, 13, 12, 10}, {9, 14, 20, 16, 11, 13, 22, 18}
452  };
453 
454  for (int j = 0; j < Nvert; j++)
455  {
456  mesh_split[3]->AddVertex(hex_v[j]);
457  }
458  for (int j = 0; j < NEsplit; j++)
459  {
460  int attribute = j + 1;
461  mesh_split[3]->AddHex(hex_e[j], attribute);
462  }
463  mesh_split[3]->FinalizeHexMesh(1, 1, true);
464 
466  gf_rst_map[3] = new GridFunction(fes_rst_map[3]);
467  const int npt = gf_rst_map[3]->Size()/dim;
468  for (int k = 0; k < dim; k++)
469  {
470  for (int j = 0; j < npt; j++)
471  {
472  (*gf_rst_map[3])(j+k*npt) = hex_v[j][k];
473  }
474  }
475  }
476  }
477 
478  NE_split_total = 0;
481  int NEsplit = 0;
482  for (int e = 0; e < mesh->GetNE(); e++)
483  {
484  const Geometry::Type gt = mesh->GetElement(e)->GetGeometryType();
485  if (gt == Geometry::TRIANGLE || gt == Geometry::PRISM)
486  {
487  NEsplit = 3;
488  }
489  else if (gt == Geometry::TETRAHEDRON)
490  {
491  NEsplit = 4;
492  }
493  else if (gt == Geometry::PYRAMID)
494  {
495  NEsplit = 8;
496  }
497  else if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
498  {
499  NEsplit = 1;
500  }
501  else
502  {
503  MFEM_ABORT("Unsupported geometry type.");
504  }
505  NE_split_total += NEsplit;
506  for (int i = 0; i < NEsplit; i++)
507  {
510  }
511  }
512 }
513 
515  IntegrationRule *irule,
516  int order)
517 {
518  H1_FECollection fec(order, dim);
519  FiniteElementSpace nodal_fes(meshin, &fec, dim);
520  meshin->SetNodalFESpace(&nodal_fes);
521  const int NEsplit = meshin->GetNE();
522 
523  const int dof_cnt = nodal_fes.GetFE(0)->GetDof(),
524  pts_cnt = NEsplit * dof_cnt;
525  Vector irlist(dim * pts_cnt);
526 
527  const TensorBasisElement *tbe =
528  dynamic_cast<const TensorBasisElement *>(nodal_fes.GetFE(0));
529  MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
530  const Array<int> &dof_map = tbe->GetDofMap();
531 
532  DenseMatrix pos(dof_cnt, dim);
533  Vector posV(pos.Data(), dof_cnt * dim);
534  Array<int> xdofs(dof_cnt * dim);
535 
536  // Create an IntegrationRule on the nodes of the reference submesh.
537  MFEM_ASSERT(irule->GetNPoints() == pts_cnt, "IntegrationRule does not have"
538  "the correct number of points.");
539  GridFunction *nodesplit = meshin->GetNodes();
540  int pt_id = 0;
541  for (int i = 0; i < NEsplit; i++)
542  {
543  nodal_fes.GetElementVDofs(i, xdofs);
544  nodesplit->GetSubVector(xdofs, posV);
545  for (int j = 0; j < dof_cnt; j++)
546  {
547  for (int d = 0; d < dim; d++)
548  {
549  irlist(pts_cnt * d + pt_id) = pos(dof_map[j], d);
550  }
551  irule->IntPoint(pt_id).x = irlist(pt_id);
552  irule->IntPoint(pt_id).y = irlist(pts_cnt + pt_id);
553  if (dim == 3)
554  {
555  irule->IntPoint(pt_id).z = irlist(2*pts_cnt + pt_id);
556  }
557  pt_id++;
558  }
559  }
560 }
561 
563  Vector &node_vals)
564 {
565  const GridFunction *nodes = gf_in;
566  const FiniteElementSpace *fes = nodes->FESpace();
567  const int NE = mesh->GetNE();
568  const int vdim = gf_in->FESpace()->GetVDim();
569 
570  IntegrationRule *ir_split_temp = NULL;
571 
572  const int dof_1D = nodes->FESpace()->GetFE(0)->GetOrder()+1;
573  const int pts_el = std::pow(dof_1D, dim);
574  const int pts_cnt = NE_split_total * pts_el;
575  node_vals.SetSize(vdim * pts_cnt);
576  node_vals = 0.0;
577 
578  int gsl_mesh_pt_index = 0;
579 
580  for (int e = 0; e < NE; e++)
581  {
582  const FiniteElement *fe = nodes->FESpace()->GetFE(e);
583  const Geometry::Type gt = fe->GetGeomType();
584  bool el_to_split = true;
585  if (gt == Geometry::TRIANGLE)
586  {
587  ir_split_temp = ir_split[0];
588  }
589  else if (gt == Geometry::TETRAHEDRON)
590  {
591  ir_split_temp = ir_split[1];
592  }
593  else if (gt == Geometry::PRISM)
594  {
595  ir_split_temp = ir_split[2];
596  }
597  else if (gt == Geometry::PYRAMID)
598  {
599  ir_split_temp = ir_split[3];
600  }
601  else if (gt == Geometry::SQUARE || gt == Geometry::CUBE)
602  {
603  el_to_split = false;
604  }
605  else
606  {
607  MFEM_ABORT("Unsupported geometry type.");
608  }
609 
610  if (el_to_split) // Triangle/Tet/Prism
611  {
612  // Fill gsl_mesh with location of split points.
613  Vector locval(vdim);
614  for (int i = 0; i < ir_split_temp->GetNPoints(); i++)
615  {
616  const IntegrationPoint &ip = ir_split_temp->IntPoint(i);
617  nodes->GetVectorValue(e, ip, locval);
618  for (int d = 0; d < vdim; d++)
619  {
620  node_vals(pts_cnt*d + gsl_mesh_pt_index) = locval(d);
621  }
622  gsl_mesh_pt_index++;
623  }
624  }
625  else // Quad/Hex
626  {
627  const int dof_cnt_split = fe->GetDof();
628 
629  const TensorBasisElement *tbe =
630  dynamic_cast<const TensorBasisElement *>(fes->GetFE(e));
631  MFEM_VERIFY(tbe != NULL, "TensorBasis FiniteElement expected.");
632  Array<int> dof_map(dof_cnt_split);
633  const Array<int> &dm = tbe->GetDofMap();
634  if (dm.Size() > 0) { dof_map = dm; }
635  else { for (int i = 0; i < dof_cnt_split; i++) { dof_map[i] = i; } }
636 
637  DenseMatrix pos(dof_cnt_split, vdim);
638  Vector posV(pos.Data(), dof_cnt_split * vdim);
639  Array<int> xdofs(dof_cnt_split * vdim);
640 
641  fes->GetElementVDofs(e, xdofs);
642  nodes->GetSubVector(xdofs, posV);
643  for (int j = 0; j < dof_cnt_split; j++)
644  {
645  for (int d = 0; d < vdim; d++)
646  {
647  node_vals(pts_cnt * d + gsl_mesh_pt_index) = pos(dof_map[j], d);
648  }
649  gsl_mesh_pt_index++;
650  }
651  }
652  }
653 }
654 
655 
657 {
660 
661  gsl_mfem_ref -= -1.; // map [-1, 1] to
662  gsl_mfem_ref *= 0.5; // [0, 1]
663 
664  int nptorig = points_cnt,
665  npt = points_cnt;
666 
667  GridFunction *gf_rst_map_temp = NULL;
668  int nptsend = 0;
669 
670  for (int index = 0; index < npt; index++)
671  {
672  if (gsl_code[index] != 2 && gsl_proc[index] != gsl_comm->id)
673  {
674  nptsend +=1;
675  }
676  }
677 
678  // Pack data to send via crystal router
679  struct gslib::array *outpt = new gslib::array;
680  struct out_pt { double r[3]; uint index, el, proc; };
681  struct out_pt *pt;
682  array_init(struct out_pt, outpt, nptsend);
683  outpt->n=nptsend;
684  pt = (struct out_pt *)outpt->ptr;
685  for (int index = 0; index < npt; index++)
686  {
687  if (gsl_code[index] == 2 || gsl_proc[index] == gsl_comm->id)
688  {
689  continue;
690  }
691  for (int d = 0; d < dim; ++d)
692  {
693  pt->r[d]= gsl_mfem_ref(index*dim + d);
694  }
695  pt->index = index;
696  pt->proc = gsl_proc[index];
697  pt->el = gsl_elem[index];
698  ++pt;
699  }
700 
701  // Transfer data to target MPI ranks
702  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
703 
704  // Map received points
705  npt = outpt->n;
706  pt = (struct out_pt *)outpt->ptr;
707  for (int index = 0; index < npt; index++)
708  {
709  IntegrationPoint ip;
710  ip.Set3(&pt->r[0]);
711  const int elem = pt->el;
712  const int mesh_elem = split_element_map[elem];
713  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(mesh_elem);
714  const Geometry::Type gt = fe->GetGeomType();
715  pt->el = mesh_elem;
716 
717  if (gt == Geometry::SQUARE || gt == Geometry::CUBE) { ++pt; continue; }
718  else if (gt == Geometry::TRIANGLE)
719  {
720  gf_rst_map_temp = gf_rst_map[0];
721  }
722  else if (gt == Geometry::TETRAHEDRON)
723  {
724  gf_rst_map_temp = gf_rst_map[1];
725  }
726  else if (gt == Geometry::PRISM)
727  {
728  gf_rst_map_temp = gf_rst_map[2];
729  }
730  else if (gt == Geometry::PYRAMID)
731  {
732  gf_rst_map_temp = gf_rst_map[3];
733  }
734 
735  int local_elem = split_element_index[elem];
736  Vector mfem_ref(dim);
737  // map to rst of macro element
738  gf_rst_map_temp->GetVectorValue(local_elem, ip, mfem_ref);
739 
740  for (int d = 0; d < dim; d++)
741  {
742  pt->r[d] = mfem_ref(d);
743  }
744  ++pt;
745  }
746 
747  // Transfer data back to source MPI rank
748  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
749  npt = outpt->n;
750 
751  // First copy mapped information for points on other procs
752  pt = (struct out_pt *)outpt->ptr;
753  for (int index = 0; index < npt; index++)
754  {
755  gsl_mfem_elem[pt->index] = pt->el;
756  for (int d = 0; d < dim; d++)
757  {
758  gsl_mfem_ref(d + pt->index*dim) = pt->r[d];
759  }
760  ++pt;
761  }
762  array_free(outpt);
763  delete outpt;
764 
765  // Now map information for points on the same proc
766  for (int index = 0; index < nptorig; index++)
767  {
768  if (gsl_code[index] != 2 && gsl_proc[index] == gsl_comm->id)
769  {
770  const int elem = gsl_elem[index];
771  const int mesh_elem = split_element_map[elem];
772  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(mesh_elem);
773  const Geometry::Type gt = fe->GetGeomType();
774  gsl_mfem_elem[index] = mesh_elem;
775  if (gt == Geometry::SQUARE || gt == Geometry::CUBE) { continue; }
776  else if (gt == Geometry::TRIANGLE)
777  {
778  gf_rst_map_temp = gf_rst_map[0];
779  }
780  else if (gt == Geometry::TETRAHEDRON)
781  {
782  gf_rst_map_temp = gf_rst_map[1];
783  }
784  else if (gt == Geometry::PRISM)
785  {
786  gf_rst_map_temp = gf_rst_map[2];
787  }
788  else if (gt == Geometry::PYRAMID)
789  {
790  gf_rst_map_temp = gf_rst_map[3];
791  }
792 
793  int local_elem = split_element_index[elem];
794  IntegrationPoint ip;
795  Vector mfem_ref(gsl_mfem_ref.GetData()+index*dim, dim);
796  ip.Set2(mfem_ref.GetData());
797  if (dim == 3) { ip.z = mfem_ref(2); }
798  gf_rst_map_temp->GetVectorValue(local_elem, ip, mfem_ref);
799  }
800  }
801 }
802 
804  Vector &field_out)
805 {
806  const int gf_order = field_in.FESpace()->GetFE(0)->GetOrder(),
807  mesh_order = mesh->GetNodalFESpace()->GetFE(0)->GetOrder();
808 
809  const FiniteElementCollection *fec_in = field_in.FESpace()->FEColl();
810  const H1_FECollection *fec_h1 = dynamic_cast<const H1_FECollection *>(fec_in);
811  const L2_FECollection *fec_l2 = dynamic_cast<const L2_FECollection *>(fec_in);
812 
813  if (fec_h1 && gf_order == mesh_order &&
814  fec_h1->GetBasisType() == BasisType::GaussLobatto &&
815  !field_in.FESpace()->IsVariableOrder())
816  {
817  InterpolateH1(field_in, field_out);
818  return;
819  }
820  else
821  {
822  InterpolateGeneral(field_in, field_out);
823  if (!fec_l2 || avgtype == AvgType::NONE) { return; }
824  }
825 
826  // For points on element borders, project the L2 GridFunction to H1 and
827  // re-interpolate.
828  if (fec_l2)
829  {
830  Array<int> indl2;
831  for (int i = 0; i < points_cnt; i++)
832  {
833  if (gsl_code[i] == 1) { indl2.Append(i); }
834  }
835  int borderPts = indl2.Size();
836 #ifdef MFEM_USE_MPI
837  MPI_Allreduce(MPI_IN_PLACE, &borderPts, 1, MPI_INT, MPI_SUM, gsl_comm->c);
838 #endif
839  if (borderPts == 0) { return; } // no points on element borders
840 
841  Vector field_out_l2(field_out.Size());
842  VectorGridFunctionCoefficient field_in_dg(&field_in);
843  int gf_order_h1 = std::max(gf_order, 1); // H1 should be at least order 1
844  H1_FECollection fec(gf_order_h1, dim);
845  const int ncomp = field_in.FESpace()->GetVDim();
846  FiniteElementSpace fes(mesh, &fec, ncomp);
847  GridFunction field_in_h1(&fes);
848 
849  if (avgtype == AvgType::ARITHMETIC)
850  {
851  field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::ARITHMETIC);
852  }
853  else if (avgtype == AvgType::HARMONIC)
854  {
855  field_in_h1.ProjectDiscCoefficient(field_in_dg, GridFunction::HARMONIC);
856  }
857  else
858  {
859  MFEM_ABORT("Invalid averaging type.");
860  }
861 
862  if (gf_order_h1 == mesh_order) // basis is GaussLobatto by default
863  {
864  InterpolateH1(field_in_h1, field_out_l2);
865  }
866  else
867  {
868  InterpolateGeneral(field_in_h1, field_out_l2);
869  }
870 
871  // Copy interpolated values for the points on element border
872  for (int j = 0; j < ncomp; j++)
873  {
874  for (int i = 0; i < indl2.Size(); i++)
875  {
876  int idx = field_in.FESpace()->GetOrdering() == Ordering::byNODES ?
877  indl2[i] + j*points_cnt:
878  indl2[i]*ncomp + j;
879  field_out(idx) = field_out_l2(idx);
880  }
881  }
882  }
883 }
884 
886  Vector &field_out)
887 {
888  FiniteElementSpace ind_fes(mesh, field_in.FESpace()->FEColl());
889  GridFunction field_in_scalar(&ind_fes);
890  Vector node_vals;
891 
892  const int ncomp = field_in.FESpace()->GetVDim(),
893  points_fld = field_in.Size() / ncomp,
895 
896  field_out.SetSize(points_cnt*ncomp);
897  field_out = default_interp_value;
898 
899  for (int i = 0; i < ncomp; i++)
900  {
901  const int dataptrin = i*points_fld,
902  dataptrout = i*points_cnt;
903  if (field_in.FESpace()->GetOrdering() == Ordering::byNODES)
904  {
905  field_in_scalar.NewDataAndSize(field_in.GetData()+dataptrin, points_fld);
906  }
907  else
908  {
909  for (int j = 0; j < points_fld; j++)
910  {
911  field_in_scalar(j) = field_in(i + j*ncomp);
912  }
913  }
914  GetNodalValues(&field_in_scalar, node_vals);
915 
916  if (dim==2)
917  {
918  findpts_eval_2(field_out.GetData()+dataptrout, sizeof(double),
919  gsl_code.GetData(), sizeof(unsigned int),
920  gsl_proc.GetData(), sizeof(unsigned int),
921  gsl_elem.GetData(), sizeof(unsigned int),
922  gsl_ref.GetData(), sizeof(double) * dim,
923  points_cnt, node_vals.GetData(), fdata2D);
924  }
925  else
926  {
927  findpts_eval_3(field_out.GetData()+dataptrout, sizeof(double),
928  gsl_code.GetData(), sizeof(unsigned int),
929  gsl_proc.GetData(), sizeof(unsigned int),
930  gsl_elem.GetData(), sizeof(unsigned int),
931  gsl_ref.GetData(), sizeof(double) * dim,
932  points_cnt, node_vals.GetData(), fdata3D);
933  }
934  }
935  if (field_in.FESpace()->GetOrdering() == Ordering::byVDIM)
936  {
937  Vector field_out_temp = field_out;
938  for (int i = 0; i < ncomp; i++)
939  {
940  for (int j = 0; j < points_cnt; j++)
941  {
942  field_out(i + j*ncomp) = field_out_temp(j + i*points_cnt);
943  }
944  }
945  }
946 }
947 
949  Vector &field_out)
950 {
951  int ncomp = field_in.VectorDim(),
952  nptorig = points_cnt,
953  npt = points_cnt;
954 
955  field_out.SetSize(points_cnt*ncomp);
956  field_out = default_interp_value;
957 
958  if (gsl_comm->np == 1) // serial
959  {
960  for (int index = 0; index < npt; index++)
961  {
962  if (gsl_code[index] == 2) { continue; }
963  IntegrationPoint ip;
965  if (dim == 3) { ip.z = gsl_mfem_ref(index*dim + 2); }
966  Vector localval(ncomp);
967  field_in.GetVectorValue(gsl_mfem_elem[index], ip, localval);
968  if (field_in.FESpace()->GetOrdering() == Ordering::byNODES)
969  {
970  for (int i = 0; i < ncomp; i++)
971  {
972  field_out(index + i*npt) = localval(i);
973  }
974  }
975  else //byVDIM
976  {
977  for (int i = 0; i < ncomp; i++)
978  {
979  field_out(index*ncomp + i) = localval(i);
980  }
981  }
982  }
983  }
984  else // parallel
985  {
986  // Determine number of points to be sent
987  int nptsend = 0;
988  for (int index = 0; index < npt; index++)
989  {
990  if (gsl_code[index] != 2) { nptsend +=1; }
991  }
992 
993  // Pack data to send via crystal router
994  struct gslib::array *outpt = new gslib::array;
995  struct out_pt { double r[3], ival; uint index, el, proc; };
996  struct out_pt *pt;
997  array_init(struct out_pt, outpt, nptsend);
998  outpt->n=nptsend;
999  pt = (struct out_pt *)outpt->ptr;
1000  for (int index = 0; index < npt; index++)
1001  {
1002  if (gsl_code[index] == 2) { continue; }
1003  for (int d = 0; d < dim; ++d) { pt->r[d]= gsl_mfem_ref(index*dim + d); }
1004  pt->index = index;
1005  pt->proc = gsl_proc[index];
1006  pt->el = gsl_mfem_elem[index];
1007  ++pt;
1008  }
1009 
1010  // Transfer data to target MPI ranks
1011  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
1012 
1013  if (ncomp == 1)
1014  {
1015  // Interpolate the grid function
1016  npt = outpt->n;
1017  pt = (struct out_pt *)outpt->ptr;
1018  for (int index = 0; index < npt; index++)
1019  {
1020  IntegrationPoint ip;
1021  ip.Set3(&pt->r[0]);
1022  pt->ival = field_in.GetValue(pt->el, ip, 1);
1023  ++pt;
1024  }
1025 
1026  // Transfer data back to source MPI rank
1027  sarray_transfer(struct out_pt, outpt, proc, 1, cr);
1028  npt = outpt->n;
1029  pt = (struct out_pt *)outpt->ptr;
1030  for (int index = 0; index < npt; index++)
1031  {
1032  field_out(pt->index) = pt->ival;
1033  ++pt;
1034  }
1035  array_free(outpt);
1036  delete outpt;
1037  }
1038  else // ncomp > 1
1039  {
1040  // Interpolate data and store in a Vector
1041  npt = outpt->n;
1042  pt = (struct out_pt *)outpt->ptr;
1043  Vector vec_int_vals(npt*ncomp);
1044  for (int index = 0; index < npt; index++)
1045  {
1046  IntegrationPoint ip;
1047  ip.Set3(&pt->r[0]);
1048  Vector localval(vec_int_vals.GetData()+index*ncomp, ncomp);
1049  field_in.GetVectorValue(pt->el, ip, localval);
1050  ++pt;
1051  }
1052 
1053  // Save index and proc data in a struct
1054  struct gslib::array *savpt = new gslib::array;
1055  struct sav_pt { uint index, proc; };
1056  struct sav_pt *spt;
1057  array_init(struct sav_pt, savpt, npt);
1058  savpt->n=npt;
1059  spt = (struct sav_pt *)savpt->ptr;
1060  pt = (struct out_pt *)outpt->ptr;
1061  for (int index = 0; index < npt; index++)
1062  {
1063  spt->index = pt->index;
1064  spt->proc = pt->proc;
1065  ++pt; ++spt;
1066  }
1067 
1068  array_free(outpt);
1069  delete outpt;
1070 
1071  // Copy data from save struct to send struct and send component wise
1072  struct gslib::array *sendpt = new gslib::array;
1073  struct send_pt { double ival; uint index, proc; };
1074  struct send_pt *sdpt;
1075  for (int j = 0; j < ncomp; j++)
1076  {
1077  array_init(struct send_pt, sendpt, npt);
1078  sendpt->n=npt;
1079  spt = (struct sav_pt *)savpt->ptr;
1080  sdpt = (struct send_pt *)sendpt->ptr;
1081  for (int index = 0; index < npt; index++)
1082  {
1083  sdpt->index = spt->index;
1084  sdpt->proc = spt->proc;
1085  sdpt->ival = vec_int_vals(j + index*ncomp);
1086  ++sdpt; ++spt;
1087  }
1088 
1089  sarray_transfer(struct send_pt, sendpt, proc, 1, cr);
1090  sdpt = (struct send_pt *)sendpt->ptr;
1091  for (int index = 0; index < static_cast<int>(sendpt->n); index++)
1092  {
1093  int idx = field_in.FESpace()->GetOrdering() == Ordering::byNODES ?
1094  sdpt->index + j*nptorig :
1095  sdpt->index*ncomp + j;
1096  field_out(idx) = sdpt->ival;
1097  ++sdpt;
1098  }
1099  array_free(sendpt);
1100  }
1101  array_free(savpt);
1102  delete sendpt;
1103  delete savpt;
1104  } // ncomp > 1
1105  } // parallel
1106 }
1107 
1108 void OversetFindPointsGSLIB::Setup(Mesh &m, const int meshid,
1109  GridFunction *gfmax,
1110  const double bb_t, const double newt_tol,
1111  const int npt_max)
1112 {
1113  MFEM_VERIFY(m.GetNodes() != NULL, "Mesh nodes are required.");
1114  MFEM_VERIFY(!(m.GetNodes()->FESpace()->IsVariableOrder()),
1115  "Variable order mesh is not currently supported.");
1116 
1117  // FreeData if OversetFindPointsGSLIB::Setup has been called already
1118  if (setupflag) { FreeData(); }
1119 
1120  crystal_init(cr, gsl_comm);
1121  mesh = &m;
1122  dim = mesh->Dimension();
1123  const FiniteElement *fe = mesh->GetNodalFESpace()->GetFE(0);
1124  unsigned dof1D = fe->GetOrder() + 1;
1125 
1126  SetupSplitMeshes();
1127  if (dim == 2)
1128  {
1129  if (ir_split[0]) { delete ir_split[0]; ir_split[0] = NULL; }
1130  ir_split[0] = new IntegrationRule(3*pow(dof1D, dim));
1132  }
1133  else if (dim == 3)
1134  {
1135  if (ir_split[1]) { delete ir_split[1]; ir_split[1] = NULL; }
1136  ir_split[1] = new IntegrationRule(4*pow(dof1D, dim));
1138 
1139  if (ir_split[2]) { delete ir_split[2]; ir_split[2] = NULL; }
1140  ir_split[2] = new IntegrationRule(3*pow(dof1D, dim));
1142 
1143  if (ir_split[3]) { delete ir_split[3]; ir_split[3] = NULL; }
1144  ir_split[3] = new IntegrationRule(8*pow(dof1D, dim));
1146  }
1147 
1149 
1150  MFEM_ASSERT(meshid>=0, " The ID should be greater than or equal to 0.");
1151 
1152  const int pts_cnt = gsl_mesh.Size()/dim,
1153  NEtot = pts_cnt/(int)pow(dof1D, dim);
1154 
1155  distfint.SetSize(pts_cnt);
1156  if (!gfmax)
1157  {
1158  distfint = 0.0;
1159  }
1160  else
1161  {
1162  GetNodalValues(gfmax, distfint);
1163  }
1164  u_meshid = (unsigned int)meshid;
1165 
1166  if (dim == 2)
1167  {
1168  unsigned nr[2] = { dof1D, dof1D };
1169  unsigned mr[2] = { 2*dof1D, 2*dof1D };
1170  double * const elx[2] = { &gsl_mesh(0), &gsl_mesh(pts_cnt) };
1171  fdata2D = findptsms_setup_2(gsl_comm, elx, nr, NEtot, mr, bb_t,
1172  pts_cnt, pts_cnt, npt_max, newt_tol,
1173  &u_meshid, &distfint(0));
1174  }
1175  else
1176  {
1177  unsigned nr[3] = { dof1D, dof1D, dof1D };
1178  unsigned mr[3] = { 2*dof1D, 2*dof1D, 2*dof1D };
1179  double * const elx[3] =
1180  { &gsl_mesh(0), &gsl_mesh(pts_cnt), &gsl_mesh(2*pts_cnt) };
1181  fdata3D = findptsms_setup_3(gsl_comm, elx, nr, NEtot, mr, bb_t,
1182  pts_cnt, pts_cnt, npt_max, newt_tol,
1183  &u_meshid, &distfint(0));
1184  }
1185  setupflag = true;
1186  overset = true;
1187 }
1188 
1190  Array<unsigned int> &point_id,
1191  int point_pos_ordering)
1192 {
1193  MFEM_VERIFY(setupflag, "Use OversetFindPointsGSLIB::Setup before "
1194  "finding points.");
1195  MFEM_VERIFY(overset, " Please setup FindPoints for overlapping grids.");
1196  points_cnt = point_pos.Size() / dim;
1197  unsigned int match = 0; // Don't find points in the mesh if point_id = mesh_id
1198 
1204 
1205  auto xvFill = [&](const double *xv_base[], unsigned xv_stride[], int dim)
1206  {
1207  for (int d = 0; d < dim; d++)
1208  {
1209  if (point_pos_ordering == Ordering::byNODES)
1210  {
1211  xv_base[d] = point_pos.GetData() + d*points_cnt;
1212  xv_stride[d] = sizeof(double);
1213  }
1214  else
1215  {
1216  xv_base[d] = point_pos.GetData() + d;
1217  xv_stride[d] = dim*sizeof(double);
1218  }
1219  }
1220  };
1221  if (dim == 2)
1222  {
1223  const double *xv_base[2];
1224  unsigned xv_stride[2];
1225  xvFill(xv_base, xv_stride, dim);
1226  findptsms_2(gsl_code.GetData(), sizeof(unsigned int),
1227  gsl_proc.GetData(), sizeof(unsigned int),
1228  gsl_elem.GetData(), sizeof(unsigned int),
1229  gsl_ref.GetData(), sizeof(double) * dim,
1230  gsl_dist.GetData(), sizeof(double),
1231  xv_base, xv_stride,
1232  point_id.GetData(), sizeof(unsigned int), &match,
1233  points_cnt, fdata2D);
1234  }
1235  else // dim == 3
1236  {
1237  const double *xv_base[3];
1238  unsigned xv_stride[3];
1239  xvFill(xv_base, xv_stride, dim);
1240  findptsms_3(gsl_code.GetData(), sizeof(unsigned int),
1241  gsl_proc.GetData(), sizeof(unsigned int),
1242  gsl_elem.GetData(), sizeof(unsigned int),
1243  gsl_ref.GetData(), sizeof(double) * dim,
1244  gsl_dist.GetData(), sizeof(double),
1245  xv_base, xv_stride,
1246  point_id.GetData(), sizeof(unsigned int), &match,
1247  points_cnt, fdata3D);
1248  }
1249 
1250  // Set the element number and reference position to 0 for points not found
1251  for (int i = 0; i < points_cnt; i++)
1252  {
1253  if (gsl_code[i] == 2 ||
1254  (gsl_code[i] == 1 && gsl_dist(i) > bdr_tol))
1255  {
1256  gsl_elem[i] = 0;
1257  for (int d = 0; d < dim; d++) { gsl_ref(i*dim + d) = -1.; }
1258  gsl_code[i] = 2;
1259  }
1260  }
1261 
1262  // Map element number for simplices, and ref_pos from [-1,1] to [0,1] for both
1263  // simplices and quads.
1265 }
1266 
1268  Array<unsigned int> &point_id,
1269  const GridFunction &field_in,
1270  Vector &field_out,
1271  int point_pos_ordering)
1272 {
1273  FindPoints(point_pos, point_id, point_pos_ordering);
1274  Interpolate(field_in, field_out);
1275 }
1276 
1277 
1278 } // namespace mfem
1279 
1280 #endif // MFEM_USE_GSLIB
Abstract class for all finite elements.
Definition: fe_base.hpp:232
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:803
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Array< unsigned int > gsl_elem
Definition: gslib.hpp:70
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5351
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:463
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:479
virtual ~FindPointsGSLIB()
Definition: gslib.cpp:67
void FindPoints(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:171
int Dimension() const
Definition: mesh.hpp:1047
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:512
Array< GridFunction * > gf_rst_map
Definition: gslib.hpp:63
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2570
Array< unsigned int > gsl_mfem_elem
Definition: gslib.hpp:70
T * GetData()
Returns the data.
Definition: array.hpp:115
unsigned int uint
Definition: gecko.hpp:204
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
virtual void SetupSplitMeshes()
Definition: gslib.cpp:296
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
double * Data() const
Returns the matrix data array.
Definition: densemat.hpp:111
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1267
const Element * GetElement(int i) const
Definition: mesh.hpp:1081
int GetBasisType() const
Definition: fe_coll.hpp:275
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:548
Array< unsigned int > gsl_proc
Definition: gslib.hpp:70
Array< Mesh * > mesh_split
Definition: gslib.hpp:59
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:452
void DeleteAll()
Delete the whole array.
Definition: array.hpp:851
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
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe_base.hpp:1188
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:319
virtual void MapRefPosAndElemIndices()
Definition: gslib.cpp:656
struct gslib::findpts_data_3 * fdata3D
Definition: gslib.hpp:66
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals)
Get GridFunction value at the points expected by GSLIB.
Definition: gslib.cpp:562
Array< IntegrationRule * > ir_split
Definition: gslib.hpp:60
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:756
virtual void FreeData()
Definition: gslib.cpp:267
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:601
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:107
double default_interp_value
Definition: gslib.hpp:73
struct gslib::crystal * cr
Definition: gslib.hpp:67
void Set2(const double x1, const double x2)
Definition: intrules.hpp:80
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5298
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:683
Array< int > split_element_index
Definition: gslib.hpp:76
Array< FiniteElementSpace * > fes_rst_map
Definition: gslib.hpp:62
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:581
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
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:70
Array< unsigned int > gsl_code
Definition: gslib.hpp:70
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:948
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:684
virtual void SetupIntegrationRuleForSplitMesh(Mesh *mesh, IntegrationRule *irule, int order)
Definition: gslib.cpp:514
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:322
struct gslib::findpts_data_2 * fdata2D
Definition: gslib.hpp:65
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:936
Class for integration point with weight.
Definition: intrules.hpp:25
int VectorDim() const
Definition: gridfunc.cpp:324
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:599
Array< int > split_element_map
Definition: gslib.hpp:75
int dim
Definition: ex24.cpp:53
void Destroy()
Destroy a vector.
Definition: vector.hpp:589
int index(int i, int j, int nx, int ny)
Definition: life.cpp:235
struct gslib::comm * gsl_comm
Definition: gslib.hpp:68
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
FiniteElementCollection * fec_map_lin
Definition: gslib.hpp:64
Vector data type.
Definition: vector.hpp:60
Vector coefficient defined by a vector GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:252
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7949
Definition: gslib.cpp:23
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1189
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:1108
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:326
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition: gslib.cpp:885
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:320