MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pumi.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 "pumi.hpp"
13 
14 #ifdef MFEM_USE_PUMI
15 #ifdef MFEM_USE_MPI
16 
17 #include "mesh_headers.hpp"
18 #include "../fem/fem.hpp"
19 #include "../general/sort_pairs.hpp"
20 #include "../general/text.hpp"
21 #include "../general/sets.hpp"
22 
23 #include <iostream>
24 #include <sstream>
25 #include <fstream>
26 #include <limits>
27 #include <cmath>
28 #include <cstring>
29 #include <ctime>
30 
31 using namespace std;
32 
33 namespace mfem
34 {
35 
36 static void ReadPumiElement(apf::MeshEntity* Ent, /* ptr to pumi entity */
37  apf::Downward Verts,
38  const int Attr, apf::Numbering* vert_num,
39  Element* el /* ptr to mfem entity being created */
40  )
41 {
42  int nv, *v;
43 
44  // Create element in MFEM
45  nv = el->GetNVertices();
46  v = el->GetVertices();
47 
48  // Fill the connectivity
49  for (int i = 0; i < nv; ++i)
50  {
51  v[i] = apf::getNumber(vert_num, Verts[i], 0, 0);
52  }
53 
54  // Assign attribute
55  el->SetAttribute(Attr);
56 }
57 
58 PumiMesh::PumiMesh(apf::Mesh2* apf_mesh, int generate_edges, int refine,
59  bool fix_orientation)
60 {
61  Load(apf_mesh, generate_edges, refine, fix_orientation);
62 }
63 
64 
65 
66 void PumiMesh::CountBoundaryEntity(apf::Mesh2* apf_mesh, const int BcDim,
67  int &NumBc)
68 {
69  apf::MeshEntity* ent;
70  apf::MeshIterator* itr = apf_mesh->begin(BcDim);
71 
72  while ((ent=apf_mesh->iterate(itr)))
73  {
74  apf::ModelEntity* mdEnt = apf_mesh->toModel(ent);
75  if (apf_mesh->getModelType(mdEnt) == BcDim)
76  {
77  NumBc++;
78  }
79  }
80  apf_mesh->end(itr);
81 
82  // Check if any boundary is detected
83  if (NumBc==0)
84  {
85  MFEM_ABORT("no boundary detected!");
86  }
87 }
88 
89 void PumiMesh::Load(apf::Mesh2* apf_mesh, int generate_edges, int refine,
90  bool fix_orientation)
91 {
92  int curved = 0, read_gf = 1;
93 
94  // Add a check on apf_mesh just in case
95  Clear();
96 
97  // First number vertices
98  apf::Field* apf_field_crd = apf_mesh->getCoordinateField();
99  apf::FieldShape* crd_shape = apf::getShape(apf_field_crd);
100  apf::Numbering* v_num_loc = apf::createNumbering(apf_mesh, "VertexNumbering",
101  crd_shape, 1);
102  // Check if it is a curved mesh
103  curved = (crd_shape->getOrder() > 1) ? 1 : 0;
104 
105  // Read mesh
106  ReadSCORECMesh(apf_mesh, v_num_loc, curved);
107 #ifdef MFEM_DEBUG
108  mfem::out << "After ReadSCORECMesh" << endl;
109 #endif
110  // at this point the following should be defined:
111  // 1) Dim
112  // 2) NumOfElements, elements
113  // 3) NumOfBdrElements, boundary
114  // 4) NumOfVertices, with allocated space in vertices
115  // 5) curved
116  // 5a) if curved == 0, vertices must be defined
117  // 5b) if curved != 0 and read_gf != 0,
118  // 'input' must point to a GridFunction
119  // 5c) if curved != 0 and read_gf == 0,
120  // vertices and Nodes must be defined
121 
122  // FinalizeTopology() will:
123  // - assume that generate_edges is true
124  // - assume that refine is false
125  // - does not check the orientation of regular and boundary elements
126  FinalizeTopology();
127 
128  if (curved && read_gf)
129  {
130  // Check it to be only Quadratic if higher order
131  Nodes = new GridFunctionPumi(this, apf_mesh, v_num_loc,
132  crd_shape->getOrder());
133  edge_vertex = NULL;
134  own_nodes = 1;
135  spaceDim = Nodes->VectorDim();
136 
137  // Set the 'vertices' from the 'Nodes'
138  for (int i = 0; i < spaceDim; i++)
139  {
140  Vector vert_val;
141  Nodes->GetNodalValues(vert_val, i+1);
142  for (int j = 0; j < NumOfVertices; j++)
143  {
144  vertices[j](i) = vert_val(j);
145  }
146  }
147  }
148 
149  // Delete numbering
150  apf::destroyNumbering(v_num_loc);
151 
152  Finalize(refine, fix_orientation);
153 }
154 
155 void PumiMesh::ReadSCORECMesh(apf::Mesh2* apf_mesh, apf::Numbering* v_num_loc,
156  const int curved)
157 {
158  // Here fill the element table from SCOREC MESH
159  // The vector of element pointers is generated with attr and connectivity
160 
161  apf::MeshIterator* itr = apf_mesh->begin(0);
162  apf::MeshEntity* ent;
163  NumOfVertices = 0;
164  while ((ent = apf_mesh->iterate(itr)))
165  {
166  // IDs start from 0
167  apf::number(v_num_loc, ent, 0, 0, NumOfVertices);
168  NumOfVertices++;
169  }
170  apf_mesh->end(itr);
171 
172  Dim = apf_mesh->getDimension();
173  NumOfElements = countOwned(apf_mesh,Dim);
174  elements.SetSize(NumOfElements);
175 
176  // Read elements from SCOREC Mesh
177  itr = apf_mesh->begin(Dim);
178  unsigned int j=0;
179  while ((ent = apf_mesh->iterate(itr)))
180  {
181  // Get vertices
182  apf::Downward verts;
183  apf_mesh->getDownward(ent,0,verts); // num_vert
184  // Get attribute Tag vs Geometry
185  int attr = 1;
186 
187  int geom_type = apf_mesh->getType(ent);
188  elements[j] = NewElement(geom_type);
189  ReadPumiElement(ent, verts, attr, v_num_loc, elements[j]);
190  j++;
191  }
192  // End iterator
193  apf_mesh->end(itr);
194 
195  // Read Boundaries from SCOREC Mesh
196  // First we need to count them
197  int BCdim = Dim - 1;
198  NumOfBdrElements = 0;
199  CountBoundaryEntity(apf_mesh, BCdim, NumOfBdrElements);
200  boundary.SetSize(NumOfBdrElements);
201  j=0;
202 
203  // Read boundary from SCOREC mesh
204  itr = apf_mesh->begin(BCdim);
205  while ((ent = apf_mesh->iterate(itr)))
206  {
207  // Check if this mesh entity is on the model boundary
208  apf::ModelEntity* mdEnt = apf_mesh->toModel(ent);
209  if (apf_mesh->getModelType(mdEnt) == BCdim)
210  {
211  apf::Downward verts;
212  apf_mesh->getDownward(ent, 0, verts);
213  int attr = 1;
214  int geom_type = apf_mesh->getType(ent);
215  boundary[j] = NewElement(geom_type);
216  ReadPumiElement(ent, verts, attr, v_num_loc, boundary[j]);
217  j++;
218  }
219  }
220  apf_mesh->end(itr);
221 
222  // Fill vertices
223  vertices.SetSize(NumOfVertices);
224 
225  if (!curved)
226  {
227  apf::MeshIterator* itr = apf_mesh->begin(0);
228  spaceDim = Dim;
229 
230  while ((ent = apf_mesh->iterate(itr)))
231  {
232  unsigned int id = apf::getNumber(v_num_loc, ent, 0, 0);
233  apf::Vector3 Crds;
234  apf_mesh->getPoint(ent,0,Crds);
235 
236  for (int ii=0; ii<spaceDim; ii++)
237  {
238  vertices[id](ii) = Crds[ii];
239  }
240  }
241  apf_mesh->end(itr);
242  }
243 }
244 
245 // ParPumiMesh implementation
246 // This function loads a parallel PUMI mesh and returns the parallel MFEM mesh
247 // corresponding to it.
248 ParPumiMesh::ParPumiMesh(MPI_Comm comm, apf::Mesh2* apf_mesh,
249  int refine, bool fix_orientation)
250 {
251  // Set the communicator for gtopo
252  gtopo.SetComm(comm);
253 
254  MyComm = comm;
255  MPI_Comm_size(MyComm, &NRanks);
256  MPI_Comm_rank(MyComm, &MyRank);
257 
258  Dim = apf_mesh->getDimension();
259  spaceDim = Dim;// mesh.spaceDim;
260 
261  apf::MeshIterator* itr;
262  apf::MeshEntity* ent;
263 
264  // Global numbering of vertices. This is necessary to build a local numbering
265  // that has the same ordering in each process.
266  apf::Numbering* vLocNum =
267  apf::numberOwnedDimension(apf_mesh, "AuxVertexNumbering", 0);
268  apf::GlobalNumbering* VertexNumbering = apf::makeGlobal(vLocNum, true);
269  apf::synchronize(VertexNumbering);
270 
271  // Take this process global vertex IDs and sort
272  Array<Pair<long,int> > thisVertIds(apf_mesh->count(0));
273  itr = apf_mesh->begin(0);
274  for (int i = 0; (ent = apf_mesh->iterate(itr)); i++)
275  {
276  long id = apf::getNumber(VertexNumbering, ent, 0, 0);
277  thisVertIds[i] = Pair<long,int>(id, i);
278  }
279  apf_mesh->end(itr);
280  apf::destroyGlobalNumbering(VertexNumbering);
281  thisVertIds.Sort();
282  // Set thisVertIds[i].one = j where j is such that thisVertIds[j].two = i.
283  // Thus, the mapping i -> thisVertIds[i].one is the inverse of the mapping
284  // j -> thisVertIds[j].two.
285  for (int j = 0; j < thisVertIds.Size(); j++)
286  {
287  const int i = thisVertIds[j].two;
288  thisVertIds[i].one = j;
289  }
290 
291  // Create local numbering that respects the global ordering
292  apf::Field* apf_field_crd = apf_mesh->getCoordinateField();
293  apf::FieldShape* crd_shape = apf::getShape(apf_field_crd);
294  // v_num_loc might already be associated the mesh. In that case
295  // there is no need to create it again.
296  v_num_loc = apf_mesh->findNumbering("LocalVertexNumbering");
297  if (!v_num_loc)
298  v_num_loc = apf::createNumbering(apf_mesh,
299  "LocalVertexNumbering",
300  crd_shape, 1);
301 
302  // Construct the numbering v_num_loc and set the coordinates of the vertices.
303  NumOfVertices = thisVertIds.Size();
304  vertices.SetSize(NumOfVertices);
305  itr = apf_mesh->begin(0);
306  for (int i = 0; (ent = apf_mesh->iterate(itr)); i++)
307  {
308  const int id = thisVertIds[i].one;
309  // Assign as local number
310  apf::number(v_num_loc, ent, 0, 0, id);
311 
312  apf::Vector3 Crds;
313  apf_mesh->getPoint(ent,0,Crds);
314 
315  for (int ii=0; ii<spaceDim; ii++)
316  {
317  vertices[id](ii) = Crds[ii]; // Assuming the IDs are ordered and from 0
318  }
319  }
320  apf_mesh->end(itr);
321  thisVertIds.DeleteAll();
322 
323  // Fill the elements
324  NumOfElements = countOwned(apf_mesh,Dim);
325  elements.SetSize(NumOfElements);
326 
327  // Read elements from SCOREC Mesh
328  itr = apf_mesh->begin(Dim);
329  for (int j = 0; (ent = apf_mesh->iterate(itr)); j++)
330  {
331  // Get vertices
332  apf::Downward verts;
333  apf_mesh->getDownward(ent,0,verts);
334 
335  // Get attribute Tag vs Geometry
336  int attr = 1;
337  int geom_type = apf_mesh->getType(ent);
338  elements[j] = NewElement(geom_type);
339  ReadPumiElement(ent, verts, attr, v_num_loc, elements[j]);
340  }
341  // End iterator
342  apf_mesh->end(itr);
343 
344  // Count number of boundaries by classification
345  int BcDim = Dim - 1;
346  itr = apf_mesh->begin(BcDim);
347  NumOfBdrElements = 0;
348  while ((ent=apf_mesh->iterate(itr)))
349  {
350  apf::ModelEntity* mdEnt = apf_mesh->toModel(ent);
351  if (apf_mesh->getModelType(mdEnt) == BcDim)
352  {
353  NumOfBdrElements++;
354  }
355  }
356  apf_mesh->end(itr);
357 
358  boundary.SetSize(NumOfBdrElements);
359  // Read boundary from SCOREC mesh
360  itr = apf_mesh->begin(BcDim);
361  for (int bdr_ctr = 0; (ent = apf_mesh->iterate(itr)); )
362  {
363  // Check if this mesh entity is on the model boundary
364  apf::ModelEntity* mdEnt = apf_mesh->toModel(ent);
365  if (apf_mesh->getModelType(mdEnt) == BcDim)
366  {
367  apf::Downward verts;
368  apf_mesh->getDownward(ent, 0, verts);
369  int attr = 1 ;
370  int geom_type = apf_mesh->getType(ent);
371  boundary[bdr_ctr] = NewElement(geom_type);
372  ReadPumiElement(ent, verts, attr, v_num_loc, boundary[bdr_ctr]);
373  bdr_ctr++;
374  }
375  }
376  apf_mesh->end(itr);
377 
378  // The next two methods are called by FinalizeTopology() called below:
379  // Mesh::SetMeshGen();
380  // Mesh::SetAttributes();
381 
382  // This is called by the default Mesh constructor
383  // Mesh::InitTables();
384 
385  this->FinalizeTopology();
386 
387  ListOfIntegerSets groups;
388  IntegerSet group;
389 
390  // The first group is the local one
391  group.Recreate(1, &MyRank);
392  groups.Insert(group);
393 
394  MFEM_ASSERT(Dim >= 3 || GetNFaces() == 0,
395  "[proc " << MyRank << "]: invalid state");
396 
397  // Determine shared faces
399  // Initially sfaces[i].one holds the global face id.
400  // Then it is replaced by the group id of the shared face.
401  if (Dim > 2)
402  {
403  // Number the faces globally and enumerate the local shared faces
404  // following the global enumeration. This way we ensure that the ordering
405  // of the shared faces within each group (of processors) is the same in
406  // each processor in the group.
407  apf::Numbering* AuxFaceNum =
408  apf::numberOwnedDimension(apf_mesh, "AuxFaceNumbering", 2);
409  apf::GlobalNumbering* GlobalFaceNum = apf::makeGlobal(AuxFaceNum, true);
410  apf::synchronize(GlobalFaceNum);
411 
412  itr = apf_mesh->begin(2);
413  while ((ent = apf_mesh->iterate(itr)))
414  {
415  if (apf_mesh->isShared(ent))
416  {
417  long id = apf::getNumber(GlobalFaceNum, ent, 0, 0);
418  sfaces.Append(Pair<long,apf::MeshEntity*>(id, ent));
419  }
420  }
421  apf_mesh->end(itr);
422  sfaces.Sort();
423  apf::destroyGlobalNumbering(GlobalFaceNum);
424 
425  // Replace the global face id in sfaces[i].one with group id.
426  for (int i = 0; i < sfaces.Size(); i++)
427  {
428  ent = sfaces[i].two;
429 
430  const int thisNumAdjs = 2;
431  int eleRanks[thisNumAdjs];
432 
433  // Get the IDs
434  apf::Parts res;
435  apf_mesh->getResidence(ent, res);
436  int kk = 0;
437  for (std::set<int>::iterator itr = res.begin();
438  itr != res.end(); ++itr)
439  {
440  eleRanks[kk++] = *itr;
441  }
442 
443  group.Recreate(2, eleRanks);
444  sfaces[i].one = groups.Insert(group) - 1;
445  }
446  }
447 
448  // Determine shared edges
450  // Initially sedges[i].one holds the global edge id.
451  // Then it is replaced by the group id of the shared edge.
452  if (Dim > 1)
453  {
454  // Number the edges globally and enumerate the local shared edges
455  // following the global enumeration. This way we ensure that the ordering
456  // of the shared edges within each group (of processors) is the same in
457  // each processor in the group.
458  apf::Numbering* AuxEdgeNum =
459  apf::numberOwnedDimension(apf_mesh, "EdgeNumbering", 1);
460  apf::GlobalNumbering* GlobalEdgeNum = apf::makeGlobal(AuxEdgeNum, true);
461  apf::synchronize(GlobalEdgeNum);
462 
463  itr = apf_mesh->begin(1);
464  while ((ent = apf_mesh->iterate(itr)))
465  {
466  if (apf_mesh->isShared(ent))
467  {
468  long id = apf::getNumber(GlobalEdgeNum, ent, 0, 0);
469  sedges.Append(Pair<long,apf::MeshEntity*>(id, ent));
470  }
471  }
472  apf_mesh->end(itr);
473  sedges.Sort();
474  apf::destroyGlobalNumbering(GlobalEdgeNum);
475 
476  // Replace the global edge id in sedges[i].one with group id.
477  Array<int> eleRanks;
478  for (int i = 0; i < sedges.Size(); i++)
479  {
480  ent = sedges[i].two;
481 
482  // Number of adjacent element
483  apf::Parts res;
484  apf_mesh->getResidence(ent, res);
485  eleRanks.SetSize(res.size());
486 
487  // Get the IDs
488  int kk = 0;
489  for (std::set<int>::iterator itr = res.begin();
490  itr != res.end(); itr++)
491  {
492  eleRanks[kk++] = *itr;
493  }
494 
495  // Generate the group
496  group.Recreate(eleRanks.Size(), eleRanks);
497  sedges[i].one = groups.Insert(group) - 1;
498  }
499  }
500 
501  // Determine shared vertices
503  // The entries sverts[i].one hold the local vertex ids.
504  Array<int> svert_group;
505  {
506  itr = apf_mesh->begin(0);
507  while ((ent = apf_mesh->iterate(itr)))
508  {
509  if (apf_mesh->isShared(ent))
510  {
511  int vtId = apf::getNumber(v_num_loc, ent, 0, 0);
512  sverts.Append(Pair<int,apf::MeshEntity*>(vtId, ent));
513  }
514  }
515  apf_mesh->end(itr);
516  sverts.Sort();
517 
518  // Determine svert_group
519  svert_group.SetSize(sverts.Size());
520  Array<int> eleRanks;
521  for (int i = 0; i < sverts.Size(); i++)
522  {
523  ent = sverts[i].two;
524 
525  // Number of adjacent element
526  apf::Parts res;
527  apf_mesh->getResidence(ent, res);
528  eleRanks.SetSize(res.size());
529 
530  // Get the IDs
531  int kk = 0;
532  for (std::set<int>::iterator itr = res.begin();
533  itr != res.end(); itr++)
534  {
535  eleRanks[kk++] = *itr;
536  }
537 
538  group.Recreate(eleRanks.Size(), eleRanks);
539  svert_group[i] = groups.Insert(group) - 1;
540  }
541  }
542 
543  // Build group_stria and group_squad.
544  // Also allocate shared_trias, shared_quads, and sface_lface.
545  group_stria.MakeI(groups.Size()-1);
546  group_squad.MakeI(groups.Size()-1);
547  for (int i = 0; i < sfaces.Size(); i++)
548  {
549  apf::Mesh::Type ftype = apf_mesh->getType(sfaces[i].two);
550  if (ftype == apf::Mesh::TRIANGLE)
551  {
552  group_stria.AddAColumnInRow(sfaces[i].one);
553  }
554  else if (ftype == apf::Mesh::QUAD)
555  {
556  group_squad.AddAColumnInRow(sfaces[i].one);
557  }
558  }
559  group_stria.MakeJ();
560  group_squad.MakeJ();
561  {
562  int nst = 0;
563  for (int i = 0; i < sfaces.Size(); i++)
564  {
565  apf::Mesh::Type ftype = apf_mesh->getType(sfaces[i].two);
566  if (ftype == apf::Mesh::TRIANGLE)
567  {
568  group_stria.AddConnection(sfaces[i].one, nst++);
569  }
570  else if (ftype == apf::Mesh::QUAD)
571  {
572  group_squad.AddConnection(sfaces[i].one, i-nst);
573  }
574  }
575  shared_trias.SetSize(nst);
576  shared_quads.SetSize(sfaces.Size()-nst);
577  sface_lface.SetSize(sfaces.Size());
578  }
579  group_stria.ShiftUpI();
580  group_squad.ShiftUpI();
581 
582  // Build group_sedge
583  group_sedge.MakeI(groups.Size()-1);
584  for (int i = 0; i < sedges.Size(); i++)
585  {
586  group_sedge.AddAColumnInRow(sedges[i].one);
587  }
588  group_sedge.MakeJ();
589  for (int i = 0; i < sedges.Size(); i++)
590  {
591  group_sedge.AddConnection(sedges[i].one, i);
592  }
593  group_sedge.ShiftUpI();
594 
595  // Build group_svert
596  group_svert.MakeI(groups.Size()-1);
597  for (int i = 0; i < svert_group.Size(); i++)
598  {
599  group_svert.AddAColumnInRow(svert_group[i]);
600  }
601  group_svert.MakeJ();
602  for (int i = 0; i < svert_group.Size(); i++)
603  {
604  group_svert.AddConnection(svert_group[i], i);
605  }
606  group_svert.ShiftUpI();
607 
608  // Build shared_trias and shared_quads. They are allocated above.
609  {
610  int nst = 0;
611  for (int i = 0; i < sfaces.Size(); i++)
612  {
613  ent = sfaces[i].two;
614 
615  apf::Downward verts;
616  apf_mesh->getDownward(ent,0,verts);
617 
618  int *v, nv = 0;
619  apf::Mesh::Type ftype = apf_mesh->getType(ent);
620  if (ftype == apf::Mesh::TRIANGLE)
621  {
622  v = shared_trias[nst++].v;
623  nv = 3;
624  }
625  else if (ftype == apf::Mesh::QUAD)
626  {
627  v = shared_quads[i-nst].v;
628  nv = 4;
629  }
630  for (int j = 0; j < nv; ++j)
631  {
632  v[j] = apf::getNumber(v_num_loc, verts[j], 0, 0);
633  }
634  }
635  }
636 
637  // Build shared_edges and allocate sedge_ledge
638  shared_edges.SetSize(sedges.Size());
639  sedge_ledge. SetSize(sedges.Size());
640  for (int i = 0; i < sedges.Size(); i++)
641  {
642  ent = sedges[i].two;
643 
644  apf::Downward verts;
645  apf_mesh->getDownward(ent, 0, verts);
646  int id1, id2;
647  id1 = apf::getNumber(v_num_loc, verts[0], 0, 0);
648  id2 = apf::getNumber(v_num_loc, verts[1], 0, 0);
649  if (id1 > id2) { swap(id1,id2); }
650 
651  shared_edges[i] = new Segment(id1, id2, 1);
652  }
653 
654  // Build svert_lvert
655  svert_lvert.SetSize(sverts.Size());
656  for (int i = 0; i < sverts.Size(); i++)
657  {
658  svert_lvert[i] = sverts[i].one;
659  }
660 
661  // Build the group communication topology
662  gtopo.Create(groups, 822);
663 
664  // Determine sedge_ledge and sface_lface
665  FinalizeParTopo();
666 
667  // Set nodes for higher order mesh
668  int curved = (crd_shape->getOrder() > 1) ? 1 : 0;
669  if (curved) // curved mesh
670  {
671  GridFunctionPumi auxNodes(this, apf_mesh, v_num_loc,
672  crd_shape->getOrder());
673  Nodes = new ParGridFunction(this, &auxNodes);
674  Nodes->Vector::Swap(auxNodes);
675  this->edge_vertex = NULL;
676  own_nodes = 1;
677  }
678 
679  Finalize(refine, fix_orientation);
680 }
681 
682 // GridFunctionPumi Implementation needed for high order meshes
683 GridFunctionPumi::GridFunctionPumi(Mesh* m, apf::Mesh2* PumiM,
684  apf::Numbering* v_num_loc,
685  const int mesh_order)
686 {
687  int spDim = m->SpaceDimension();
688  // Note: default BasisType for 'fec' is GaussLobatto.
689  fec = new H1_FECollection(mesh_order, m->Dimension());
690  int ordering = Ordering::byVDIM; // x1y1z1/x2y2z2/...
691  fes = new FiniteElementSpace(m, fec, spDim, ordering);
692  int data_size = fes->GetVSize();
693 
694  // Read PUMI mesh data
695  this->SetSize(data_size);
696  double* PumiData = this->GetData();
697 
698  apf::MeshEntity* ent;
699  apf::MeshIterator* itr;
700 
701  // Assume all element type are the same i.e. tetrahedral
702  const FiniteElement* H1_elem = fes->GetFE(0);
703  const IntegrationRule &All_nodes = H1_elem->GetNodes();
704  int nnodes = All_nodes.Size();
705 
706  // Loop over elements
707  apf::Field* crd_field = PumiM->getCoordinateField();
708 
709  int nc = apf::countComponents(crd_field);
710  int iel = 0;
711  itr = PumiM->begin(m->Dimension());
712  while ((ent = PumiM->iterate(itr)))
713  {
714  Array<int> vdofs;
715  fes->GetElementVDofs(iel, vdofs);
716 
717  // Create PUMI element to interpolate
718  apf::MeshElement* mE = apf::createMeshElement(PumiM, ent);
719  apf::Element* elem = apf::createElement(crd_field, mE);
720 
721  // Vertices are already interpolated
722  for (int ip = 0; ip < nnodes; ip++)
723  {
724  // Take parametric coordinates of the node
725  apf::Vector3 param;
726  param[0] = All_nodes.IntPoint(ip).x;
727  param[1] = All_nodes.IntPoint(ip).y;
728  param[2] = All_nodes.IntPoint(ip).z;
729 
730  // Compute the interpolating coordinates
731  apf::DynamicVector phCrd(nc);
732  apf::getComponents(elem, param, &phCrd[0]);
733 
734  // Fill the nodes list
735  for (int kk = 0; kk < spDim; ++kk)
736  {
737  int dof_ctr = ip + kk * nnodes;
738  PumiData[vdofs[dof_ctr]] = phCrd[kk];
739  }
740  }
741  iel++;
742  apf::destroyElement(elem);
743  apf::destroyMeshElement(mE);
744  }
745  PumiM->end(itr);
746 
747  sequence = 0;
748 }
749 
750 // Copy the adapted mesh to the original mesh and increase the sequence to be
751 // able to Call Update() methods of FESpace, Linear and Bilinear forms.
752 void ParPumiMesh::UpdateMesh(const ParMesh* AdaptedpMesh)
753 {
754  // Destroy the ParMesh data fields.
755  delete pncmesh;
756  pncmesh = NULL;
757 
758  DeleteFaceNbrData();
759 
760  for (int i = 0; i < shared_edges.Size(); i++)
761  {
762  FreeElement(shared_edges[i]);
763  }
764  shared_quads.DeleteAll();
765  shared_trias.DeleteAll();
766  shared_edges.DeleteAll();
767  group_svert.Clear();
768  group_sedge.Clear();
769  group_stria.Clear();
770  group_squad.Clear();
771  svert_lvert.DeleteAll();
772  sedge_ledge.DeleteAll();
773  sface_lface.DeleteAll();
774 
775  // Destroy the Mesh data fields.
776  Destroy();
777 
778  // Assuming Dim, spaceDim, geom type is unchanged
779  MFEM_ASSERT(Dim == AdaptedpMesh->Dim, "");
780  MFEM_ASSERT(spaceDim == AdaptedpMesh->spaceDim, "");
781  MFEM_ASSERT(meshgen == AdaptedpMesh->meshgen, "");
782 
783  NumOfVertices = AdaptedpMesh->GetNV();
784  NumOfElements = AdaptedpMesh->GetNE();
785  NumOfBdrElements = AdaptedpMesh->GetNBE();
786  NumOfEdges = AdaptedpMesh->GetNEdges();
787  NumOfFaces = AdaptedpMesh->GetNFaces();
788 
789  meshgen = AdaptedpMesh->meshgen;
790 
791  // Sequence is increased by one to trigger update in FEspace etc.
792  sequence++;
793  last_operation = Mesh::NONE;
794 
795  // Duplicate the elements
796  elements.SetSize(NumOfElements);
797  for (int i = 0; i < NumOfElements; i++)
798  {
799  elements[i] = AdaptedpMesh->GetElement(i)->Duplicate(this);
800  }
801 
802  // Copy the vertices
803  AdaptedpMesh->vertices.Copy(vertices);
804 
805  // Duplicate the boundary
806  boundary.SetSize(NumOfBdrElements);
807  for (int i = 0; i < NumOfBdrElements; i++)
808  {
809  boundary[i] = AdaptedpMesh->GetBdrElement(i)->Duplicate(this);
810  }
811 
812  // Copy the element-to-face Table, el_to_face
813  el_to_face = (AdaptedpMesh->el_to_face) ?
814  new Table(*(AdaptedpMesh->el_to_face)) : NULL;
815 
816  // Copy the boundary-to-face Array, be_to_face.
817  AdaptedpMesh->be_to_face.Copy(be_to_face);
818 
819  // Copy the element-to-edge Table, el_to_edge
820  el_to_edge = (AdaptedpMesh->el_to_edge) ?
821  new Table(*(AdaptedpMesh->el_to_edge)) : NULL;
822 
823  // Copy the boudary-to-edge Table, bel_to_edge (3D)
824  bel_to_edge = (AdaptedpMesh->bel_to_edge) ?
825  new Table(*(AdaptedpMesh->bel_to_edge)) : NULL;
826 
827  // Copy the boudary-to-edge Array, be_to_edge (2D)
828  AdaptedpMesh->be_to_edge.Copy(be_to_edge);
829 
830  // Duplicate the faces and faces_info.
831  faces.SetSize(AdaptedpMesh->faces.Size());
832  for (int i = 0; i < faces.Size(); i++)
833  {
834  Element *face = AdaptedpMesh->faces[i]; // in 1D the faces are NULL
835  faces[i] = (face) ? face->Duplicate(this) : NULL;
836  }
837  AdaptedpMesh->faces_info.Copy(faces_info);
838 
839  // Do NOT copy the element-to-element Table, el_to_el
840  el_to_el = NULL;
841 
842  // Do NOT copy the face-to-edge Table, face_edge
843  face_edge = NULL;
844 
845  // Copy the edge-to-vertex Table, edge_vertex
846  edge_vertex = (AdaptedpMesh->edge_vertex) ?
847  new Table(*(AdaptedpMesh->edge_vertex)) : NULL;
848 
849  // Copy the attributes and bdr_attributes
850  AdaptedpMesh->attributes.Copy(attributes);
851  AdaptedpMesh->bdr_attributes.Copy(bdr_attributes);
852 
853  // PUMI meshes cannot use NURBS meshes.
854  MFEM_VERIFY(AdaptedpMesh->NURBSext == NULL,
855  "invalid adapted mesh: it is a NURBS mesh");
856  NURBSext = NULL;
857 
858  // PUMI meshes cannot use NCMesh/ParNCMesh.
859  MFEM_VERIFY(AdaptedpMesh->ncmesh == NULL && AdaptedpMesh->pncmesh == NULL,
860  "invalid adapted mesh: it is a non-conforming mesh");
861  ncmesh = NULL;
862  pncmesh = NULL;
863 
864  // Parallel Implications
865  AdaptedpMesh->group_svert.Copy(group_svert);
866  AdaptedpMesh->group_sedge.Copy(group_sedge);
867  group_stria = AdaptedpMesh->group_stria;
868  group_squad = AdaptedpMesh->group_squad;
869  AdaptedpMesh->gtopo.Copy(gtopo);
870 
871  MyComm = AdaptedpMesh->MyComm;
872  NRanks = AdaptedpMesh->NRanks;
873  MyRank = AdaptedpMesh->MyRank;
874 
875  // Duplicate the shared_edges
876  shared_edges.SetSize(AdaptedpMesh->shared_edges.Size());
877  for (int i = 0; i < shared_edges.Size(); i++)
878  {
879  shared_edges[i] = AdaptedpMesh->shared_edges[i]->Duplicate(this);
880  }
881 
882  // Duplicate the shared_trias and shared_quads
883  shared_trias = AdaptedpMesh->shared_trias;
884  shared_quads = AdaptedpMesh->shared_quads;
885 
886  // Copy the shared-to-local index Arrays
887  AdaptedpMesh->svert_lvert.Copy(svert_lvert);
888  AdaptedpMesh->sedge_ledge.Copy(sedge_ledge);
889  AdaptedpMesh->sface_lface.Copy(sface_lface);
890 
891  // Do not copy face-neighbor data (can be generated if needed)
892  have_face_nbr_data = false;
893 
894  // Copy the Nodes as a ParGridFunction, including the FiniteElementCollection
895  // and the FiniteElementSpace (as a ParFiniteElementSpace)
896  if (AdaptedpMesh->Nodes)
897  {
898  FiniteElementSpace *fes = AdaptedpMesh->Nodes->FESpace();
899  const FiniteElementCollection *fec = fes->FEColl();
900  FiniteElementCollection *fec_copy =
901  FiniteElementCollection::New(fec->Name());
902  ParFiniteElementSpace *pfes_copy =
903  new ParFiniteElementSpace(this, fec_copy, fes->GetVDim(),
904  fes->GetOrdering());
905  Nodes = new ParGridFunction(pfes_copy);
906  Nodes->MakeOwner(fec_copy);
907  *Nodes = *(AdaptedpMesh->Nodes);
908  own_nodes = 1;
909  }
910 }
911 
912 int ParPumiMesh::RotationPUMItoMFEM(apf::Mesh2* apf_mesh,
913  apf::MeshEntity* tet,
914  int elemId)
915 {
916  MFEM_ASSERT(apf_mesh->getType(tet) == apf::Mesh::TET, "");
917  // get downward vertices of PUMI element
918  apf::Downward vs;
919  int nv = apf_mesh->getDownward(tet,0,vs);
920  int pumi_vid[12];
921  for (int i = 0; i < nv; i++)
922  {
923  pumi_vid[i] = apf::getNumber(v_num_loc, vs[i], 0, 0);
924  }
925 
926  // get downward vertices of MFEM element
927  mfem::Array<int> mfem_vid;
928  this->GetElementVertices(elemId, mfem_vid);
929 
930  // get rotated indices of PUMI element
931  int pumi_vid_rot[12];
932  for (int i = 0; i < nv; i++)
933  {
934  pumi_vid_rot[i] = mfem_vid.Find(pumi_vid[i]);
935  }
936  apf::Downward vs_rot;
937  for (int i = 0; i < nv; i++)
938  {
939  vs_rot[i] = vs[pumi_vid_rot[i]];
940  }
941 
942  return ma::findTetRotation(apf_mesh, tet, vs_rot);
943 }
944 
945 // Convert parent coordinate form a PUMI tet to an MFEM tet
946 IntegrationRule ParPumiMesh::ParentXisPUMItoMFEM(apf::Mesh2* apf_mesh,
947  apf::MeshEntity* tet,
948  int elemId,
949  apf::NewArray<apf::Vector3>& pumi_xi,
950  bool checkOrientation)
951 {
952  int num_nodes = pumi_xi.size();
953  IntegrationRule mfem_xi(num_nodes);
954  int rotation = checkOrientation ? RotationPUMItoMFEM(apf_mesh, tet, elemId):0;
955  for (int i = 0; i < num_nodes; i++)
956  {
957  // for non zero "rotation", rotate the xi
958  if (rotation)
959  {
960  ma::unrotateTetXi(pumi_xi[i], rotation);
961  }
962  IntegrationPoint& ip = mfem_xi.IntPoint(i);
963  double tmp_xi[3];
964  pumi_xi[i].toArray(tmp_xi);
965  ip.Set(tmp_xi,3);
966  }
967  return mfem_xi;
968 }
969 
970 // Convert parent coordinate from MFEM tet to PUMI tet
971 void ParPumiMesh::ParentXisMFEMtoPUMI(apf::Mesh2* apf_mesh,
972  int elemId,
973  apf::MeshEntity* tet,
974  const IntegrationRule& mfem_xi,
975  apf::NewArray<apf::Vector3>& pumi_xi,
976  bool checkOrientation)
977 {
978  int num_nodes = mfem_xi.Size();
979  if (!pumi_xi.allocated())
980  {
981  pumi_xi.allocate(num_nodes);
982  }
983  else
984  {
985  pumi_xi.resize(num_nodes);
986  }
987 
988  int rotation = checkOrientation ? RotationPUMItoMFEM(apf_mesh, tet, elemId):0;
989  for (int i = 0; i < num_nodes; i++)
990  {
991  IntegrationPoint ip = mfem_xi.IntPoint(i);
992  pumi_xi[i] = apf::Vector3(ip.x, ip.y, ip.z);
993 
994  // for non zero "rotation", un-rotate the xi
995  if (rotation)
996  {
997  ma::rotateTetXi(pumi_xi[i], rotation);
998  }
999  }
1000 }
1001 
1002 
1003 // Transfer a mixed vector-scalar field (i.e. velocity,pressure) and the
1004 // magnitude of the vector field to use for mesh adaptation.
1005 void ParPumiMesh::FieldMFEMtoPUMI(apf::Mesh2* apf_mesh,
1006  ParGridFunction* grid_vel,
1007  ParGridFunction* grid_pr,
1008  apf::Field* vel_field,
1009  apf::Field* pr_field,
1010  apf::Field* vel_mag_field)
1011 {
1012  apf::FieldShape* field_shape = getShape(vel_field);
1013  int dim = apf_mesh->getDimension();
1014 
1015  apf::MeshEntity* ent;
1016  apf::MeshIterator* itr = apf_mesh->begin(dim);
1017  int iel = 0;
1018  while ((ent = apf_mesh->iterate(itr)))
1019  {
1020  apf::NewArray<apf::Vector3> pumi_nodes;
1021  apf::getElementNodeXis(field_shape, apf_mesh, ent, pumi_nodes);
1022  IntegrationRule mfem_nodes = ParentXisPUMItoMFEM(
1023  apf_mesh, ent, iel, pumi_nodes, true);
1024  // Get the solution
1025  ElementTransformation* eltr = this->GetElementTransformation(iel);
1026  DenseMatrix vel;
1027  grid_vel->GetVectorValues(*eltr, mfem_nodes, vel);
1028  Vector pr;
1029  grid_pr->GetValues(iel, mfem_nodes, pr, 1);
1030 
1031  int non = 0;
1032  for (int d = 0; d <= dim; d++)
1033  {
1034  if (!field_shape->hasNodesIn(d)) { continue; }
1035  apf::Downward a;
1036  int na = apf_mesh->getDownward(ent,d,a);
1037  for (int i = 0; i < na; i++)
1038  {
1039  int type = apf_mesh->getType(a[i]);
1040  int nan = field_shape->countNodesOn(type);
1041  for (int n = 0; n < nan; n++)
1042  {
1043  apf::Vector3 v(vel.GetColumn(non));
1044  apf::setVector(vel_field, a[i], n, v);
1045  apf::setScalar(pr_field, a[i], n, pr[non]);
1046  apf::setScalar(vel_mag_field, a[i], n, v.getLength());
1047  non++;
1048  }
1049  }
1050  }
1051  iel++;
1052  }
1053  apf_mesh->end(itr);
1054 }
1055 
1056 // Transfer a scalar field its magnitude to use for mesh adaptation.
1057 void ParPumiMesh::FieldMFEMtoPUMI(apf::Mesh2* apf_mesh,
1058  ParGridFunction* grid_pr,
1059  apf::Field* pr_field,
1060  apf::Field* pr_mag_field)
1061 {
1062  apf::FieldShape* field_shape = getShape(pr_field);
1063  int dim = apf_mesh->getDimension();
1064 
1065  apf::MeshEntity* ent;
1066  apf::MeshIterator* itr = apf_mesh->begin(dim);
1067  int iel = 0;
1068  while ((ent = apf_mesh->iterate(itr)))
1069  {
1070  apf::NewArray<apf::Vector3> pumi_nodes;
1071  apf::getElementNodeXis(field_shape, apf_mesh, ent, pumi_nodes);
1072  IntegrationRule mfem_nodes = ParentXisPUMItoMFEM(
1073  apf_mesh, ent, iel, pumi_nodes, true);
1074  // Get the solution
1075  Vector vals;
1076  grid_pr->GetValues(iel, mfem_nodes, vals, 1);
1077 
1078  int non = 0;
1079  for (int d = 0; d <= dim; d++)
1080  {
1081  if (!field_shape->hasNodesIn(d)) { continue; }
1082  apf::Downward a;
1083  int na = apf_mesh->getDownward(ent,d,a);
1084  for (int i = 0; i < na; i++)
1085  {
1086  int type = apf_mesh->getType(a[i]);
1087  int nan = field_shape->countNodesOn(type);
1088  for (int n = 0; n < nan; n++)
1089  {
1090  double pr = vals[non];
1091  double pr_mag = pr >= 0 ? pr : -pr;
1092  apf::setScalar(pr_field, a[i], n, pr);
1093  apf::setScalar(pr_mag_field, a[i], n, pr_mag);
1094  non++;
1095  }
1096  }
1097  }
1098  iel++;
1099  }
1100  apf_mesh->end(itr);
1101 }
1102 
1103 // Transfer a vector field and the magnitude of the vector field to use for mesh
1104 // adaptation
1105 void ParPumiMesh::VectorFieldMFEMtoPUMI(apf::Mesh2* apf_mesh,
1106  ParGridFunction* grid_vel,
1107  apf::Field* vel_field,
1108  apf::Field* vel_mag_field)
1109 {
1110  apf::FieldShape* field_shape = getShape(vel_field);
1111  int dim = apf_mesh->getDimension();
1112 
1113  apf::MeshEntity* ent;
1114  apf::MeshIterator* itr = apf_mesh->begin(dim);
1115  int iel = 0;
1116  while ((ent = apf_mesh->iterate(itr)))
1117  {
1118  apf::NewArray<apf::Vector3> pumi_nodes;
1119  apf::getElementNodeXis(field_shape, apf_mesh, ent, pumi_nodes);
1120  IntegrationRule mfem_nodes = ParentXisPUMItoMFEM(
1121  apf_mesh, ent, iel, pumi_nodes, true);
1122  // Get the solution
1123  ElementTransformation* eltr = this->GetElementTransformation(iel);
1124  DenseMatrix vel;
1125  grid_vel->GetVectorValues(*eltr, mfem_nodes, vel);
1126 
1127  int non = 0;
1128  for (int d = 0; d <= dim; d++)
1129  {
1130  if (!field_shape->hasNodesIn(d)) { continue; }
1131  apf::Downward a;
1132  int na = apf_mesh->getDownward(ent,d,a);
1133  for (int i = 0; i < na; i++)
1134  {
1135  int type = apf_mesh->getType(a[i]);
1136  int nan = field_shape->countNodesOn(type);
1137  for (int n = 0; n < nan; n++)
1138  {
1139  apf::Vector3 v(vel.GetColumn(non));
1140  apf::setScalar(vel_mag_field, a[i], n, v.getLength());
1141  apf::setVector(vel_field, a[i], n, v);
1142  non++;
1143  }
1144  }
1145  }
1146  iel++;
1147  }
1148  apf_mesh->end(itr);
1149 }
1150 
1151 void ParPumiMesh::NedelecFieldMFEMtoPUMI(apf::Mesh2* apf_mesh,
1152  ParGridFunction* gf,
1153  apf::Field* nedelec_field)
1154 {
1155  apf::FieldShape* nedelecFieldShape = nedelec_field->getShape();
1156  int dim = apf_mesh->getDimension();
1157 
1158  // loop over all elements
1159  size_t elemNo = 0;
1160  apf::MeshEntity* ent;
1161  apf::MeshIterator* it = apf_mesh->begin(dim);
1162  while ( (ent = apf_mesh->iterate(it)) )
1163  {
1164  // get all the pumi nodes and rotate them
1165  apf::NewArray<apf::Vector3> pumi_nodes;
1166  apf::getElementNodeXis(nedelecFieldShape, apf_mesh, ent, pumi_nodes);
1167  IntegrationRule mfem_nodes = ParentXisPUMItoMFEM(
1168  apf_mesh, ent, elemNo, pumi_nodes, true);
1169  // evaluate the vector field on the mfem nodes
1170  ElementTransformation* eltr = this->GetElementTransformation(elemNo);
1171  DenseMatrix mfem_field_vals;
1172  gf->GetVectorValues(*eltr, mfem_nodes, mfem_field_vals);
1173 
1174  // compute and store dofs on ND field
1175  int non = 0;
1176  for (int d = 0; d <= dim; d++)
1177  {
1178  if (!nedelecFieldShape->hasNodesIn(d)) { continue; }
1179  apf::Downward a;
1180  int na = apf_mesh->getDownward(ent,d,a);
1181  for (int i = 0; i < na; i++)
1182  {
1183  int type = apf_mesh->getType(a[i]);
1184  int nan = nedelecFieldShape->countNodesOn(type);
1185  apf::MeshElement* me = apf::createMeshElement(apf_mesh, a[i]);
1186  for (int n = 0; n < nan; n++)
1187  {
1188  apf::Vector3 xi, tangent;
1189  nedelecFieldShape->getNodeXi(type, n, xi);
1190  nedelecFieldShape->getNodeTangent(type, n, tangent);
1191  apf::Vector3 pumi_field_vector(mfem_field_vals.GetColumn(non));
1192  apf::Matrix3x3 J;
1193  apf::getJacobian(me, xi, J);
1194  double dof = (J * pumi_field_vector) * tangent;
1195  apf::setScalar(nedelec_field, a[i], n, dof);
1196  non++;
1197  }
1198  apf::destroyMeshElement(me);
1199  }
1200  }
1201  elemNo++;
1202  }
1203  apf_mesh->end(it); // end loop over all elements
1204 }
1205 
1206 void ParPumiMesh::FieldPUMItoMFEM(apf::Mesh2* apf_mesh,
1207  apf::Field* field,
1208  ParGridFunction* grid)
1209 {
1210  int nc = apf::countComponents(field);
1211  ParFiniteElementSpace* fes = grid->ParFESpace();
1212  ParMesh* pmesh = fes->GetParMesh();
1213 
1214  int dim = apf_mesh->getDimension();
1215 
1216  apf::MeshIterator* it = apf_mesh->begin(dim);
1217  for (int i = 0; i < pmesh->GetNE(); i++)
1218  {
1219  const FiniteElement* mfem_elem = fes->GetFE(i);
1220  const IntegrationRule &mfem_xi = mfem_elem->GetNodes();
1221  int non = mfem_xi.Size();
1222  apf::MeshEntity* ent = apf_mesh->iterate(it);
1223  apf::NewArray<apf::Vector3> pumi_xi(non);
1224  ParentXisMFEMtoPUMI(apf_mesh,
1225  i,
1226  ent,
1227  mfem_xi,
1228  pumi_xi,
1229  true);
1230  Array<int> vdofs;
1231  fes->GetElementVDofs(i, vdofs);
1232  apf::MeshElement* me = apf::createMeshElement(apf_mesh, ent);
1233  apf::Element* el = apf::createElement(field, me);
1234  for (int j = 0; j < non; j++)
1235  {
1236  apf::DynamicVector values(nc);
1237  apf::getComponents(el, pumi_xi[j], &values[0]);
1238  // Fill the nodes list
1239  for (int c = 0; c < nc; c++)
1240  {
1241  int dof_loc = j + c * non;
1242  (grid->GetData())[vdofs[dof_loc]] = values[c];
1243  }
1244  }
1245  apf::destroyElement(el);
1246  apf::destroyMeshElement(me);
1247  }
1248  apf_mesh->end(it);
1249 }
1250 
1251 }
1252 
1253 #endif // MFEM_USE_MPI
1254 #endif // MFEM_USE_SCOREC
Abstract class for all finite elements.
Definition: fe.hpp:235
void Set(const double *p, const int dim)
Definition: intrules.hpp:37
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:412
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void Recreate(const int n, const int *p)
Create an integer set from C-array &#39;p&#39; of &#39;n&#39; integers. Overwrites any existing set data...
Definition: sets.cpp:59
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:633
T * end()
STL-like end. Returns pointer after the last element of the array.
Definition: array.hpp:278
int NRanks
Definition: pmesh.hpp:39
virtual Element * Duplicate(Mesh *m) const =0
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:740
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:173
int size
Size of the array.
Definition: array.hpp:50
Array< int > sface_lface
Definition: pmesh.hpp:79
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:831
Class for PUMI grid functions.
Definition: pumi.hpp:153
int VectorDim() const
Definition: gridfunc.cpp:300
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
GridFunction * Nodes
Definition: mesh.hpp:176
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
Abstract parallel finite element space.
Definition: pfespace.hpp:28
Array< Vert3 > shared_trias
Definition: pmesh.hpp:66
constexpr Element::Type QUAD
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
Array< Element * > faces
Definition: mesh.hpp:92
Array< int > sedge_ledge
Definition: pmesh.hpp:77
Table group_stria
Definition: pmesh.hpp:72
Table * el_to_face
Definition: mesh.hpp:157
ParNCMesh * pncmesh
Definition: pmesh.hpp:260
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
MPI_Comm MyComm
Definition: pmesh.hpp:38
Array< Element * > shared_edges
Definition: pmesh.hpp:62
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:726
A pair of objects.
Definition: sort_pairs.hpp:23
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe.hpp:382
int Insert(IntegerSet &s)
Check to see if set &#39;s&#39; is in the list. If not append it to the end of the list. Returns the index of...
Definition: sets.cpp:82
T * begin()
STL-like begin. Returns pointer to the first element of the array.
Definition: array.hpp:275
const Element * GetElement(int i) const
Definition: mesh.hpp:819
void Sort()
Sorts the array in ascending order. This requires operator&lt; to be defined for T.
Definition: array.hpp:234
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:455
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
int Dimension() const
Definition: mesh.hpp:788
Table * el_to_edge
Definition: mesh.hpp:156
int Size()
Return the number of integer sets in the list.
Definition: sets.hpp:67
void GetColumn(int c, Vector &col) const
Definition: densemat.cpp:1289
Array< Vert4 > shared_quads
Definition: pmesh.hpp:67
void Copy(GroupTopology &copy) const
Copy the internal data to the external &#39;copy&#39;.
int SpaceDimension() const
Definition: mesh.hpp:789
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
int Find(const T &el) const
Return the first index where &#39;el&#39; is found; return -1 if not found.
Definition: array.hpp:785
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
Array< Vertex > vertices
Definition: mesh.hpp:90
int meshgen
Definition: mesh.hpp:77
Table * bel_to_edge
Definition: mesh.hpp:160
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Definition: mesh.hpp:734
Array< int > be_to_edge
Definition: mesh.hpp:159
virtual const char * Name() const
Definition: fe_coll.hpp:61
Class for integration point with weight.
Definition: intrules.hpp:25
Table group_sedge
Definition: pmesh.hpp:71
A set of integers.
Definition: sets.hpp:23
int dim
Definition: ex24.cpp:53
Table group_svert
Shared objects in each group.
Definition: pmesh.hpp:70
Table * edge_vertex
Definition: mesh.hpp:163
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:1798
Array< FaceInfo > faces_info
Definition: mesh.hpp:153
NCMesh * ncmesh
Optional non-conforming mesh extension.
Definition: mesh.hpp:204
int Dim
Definition: mesh.hpp:66
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:743
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:746
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Array< int > svert_lvert
Shared to local index mapping.
Definition: pmesh.hpp:76
int MyRank
Definition: pmesh.hpp:39
int spaceDim
Definition: mesh.hpp:67
Class for parallel grid function.
Definition: pgridfunc.hpp:32
List of integer sets.
Definition: sets.hpp:59
Array< int > be_to_face
Definition: mesh.hpp:161
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
GroupTopology gtopo
Definition: pmesh.hpp:247
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:28
Data type line segment element.
Definition: segment.hpp:22
Table group_squad
Definition: pmesh.hpp:73
void Copy(Table &copy) const
Definition: table.cpp:390
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:199
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:823
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108