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