MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pmesh.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "mesh_headers.hpp"
17 #include "../fem/fem.hpp"
18 #include "../general/sets.hpp"
19 #include "../general/sort_pairs.hpp"
20 #include <iostream>
21 using namespace std;
22 
23 namespace mfem
24 {
25 
26 ParMesh::ParMesh(const ParMesh &pmesh, bool copy_nodes)
27  : Mesh(pmesh, false),
28  group_svert(pmesh.group_svert),
29  group_sedge(pmesh.group_sedge),
30  group_sface(pmesh.group_sface),
31  gtopo(pmesh.gtopo)
32 {
33  MyComm = pmesh.MyComm;
34  NRanks = pmesh.NRanks;
35  MyRank = pmesh.MyRank;
36 
37  // Duplicate the shared_edges
38  shared_edges.SetSize(pmesh.shared_edges.Size());
39  for (int i = 0; i < shared_edges.Size(); i++)
40  {
41  shared_edges[i] = pmesh.shared_edges[i]->Duplicate(this);
42  }
43 
44  // Duplicate the shared_faces
45  shared_faces.SetSize(pmesh.shared_faces.Size());
46  for (int i = 0; i < shared_faces.Size(); i++)
47  {
48  shared_faces[i] = pmesh.shared_faces[i]->Duplicate(this);
49  }
50 
51  // Copy the shared-to-local index Arrays
55 
56  // Do not copy face-neighbor data (can be generated if needed)
57  have_face_nbr_data = false;
58 
59  MFEM_VERIFY(pmesh.pncmesh == NULL,
60  "copying non-conforming meshes is not implemented");
61  pncmesh = NULL;
62 
63  // Copy the Nodes as a ParGridFunction, including the FiniteElementCollection
64  // and the FiniteElementSpace (as a ParFiniteElementSpace)
65  if (pmesh.Nodes && copy_nodes)
66  {
67  FiniteElementSpace *fes = pmesh.Nodes->FESpace();
68  const FiniteElementCollection *fec = fes->FEColl();
69  FiniteElementCollection *fec_copy =
71  ParFiniteElementSpace *pfes_copy =
72  new ParFiniteElementSpace(this, fec_copy, fes->GetVDim(),
73  fes->GetOrdering());
74  Nodes = new ParGridFunction(pfes_copy);
75  Nodes->MakeOwner(fec_copy);
76  *Nodes = *pmesh.Nodes;
77  own_nodes = 1;
78  }
79 }
80 
81 ParMesh::ParMesh(MPI_Comm comm, Mesh &mesh, int *partitioning_,
82  int part_method)
83  : gtopo(comm)
84 {
85  int i, j;
86  int *partitioning;
87  Array<bool> activeBdrElem;
88 
89  MyComm = comm;
90  MPI_Comm_size(MyComm, &NRanks);
91  MPI_Comm_rank(MyComm, &MyRank);
92 
93  if (mesh.Nonconforming())
94  {
95  pncmesh = new ParNCMesh(comm, *mesh.ncmesh);
96 
97  // save the element partitioning before Prune()
98  int* partition = new int[mesh.GetNE()];
99  for (int i = 0; i < mesh.GetNE(); i++)
100  {
101  partition[i] = pncmesh->InitialPartition(i);
102  }
103 
104  pncmesh->Prune();
105 
106  Mesh::Init();
109  pncmesh->OnMeshUpdated(this);
110 
111  ncmesh = pncmesh;
112  meshgen = mesh.MeshGenerator();
113 
114  mesh.attributes.Copy(attributes);
116 
118 
119  if (mesh.GetNodes())
120  {
121  Nodes = new ParGridFunction(this, mesh.GetNodes(), partition);
122  own_nodes = 1;
123  }
124  delete [] partition;
125 
126  have_face_nbr_data = false;
127  return;
128  }
129 
130  Dim = mesh.Dim;
131  spaceDim = mesh.spaceDim;
132 
133  BaseGeom = mesh.BaseGeom;
134  BaseBdrGeom = mesh.BaseBdrGeom;
135 
136  ncmesh = pncmesh = NULL;
137 
138  if (partitioning_)
139  {
140  partitioning = partitioning_;
141  }
142  else
143  {
144  partitioning = mesh.GeneratePartitioning(NRanks, part_method);
145  }
146 
147  // re-enumerate the partitions to better map to actual processor
148  // interconnect topology !?
149 
150  Array<int> vert;
151  Array<int> vert_global_local(mesh.GetNV());
152  int vert_counter, element_counter, bdrelem_counter;
153 
154  // build vert_global_local
155  vert_global_local = -1;
156 
157  element_counter = 0;
158  vert_counter = 0;
159  for (i = 0; i < mesh.GetNE(); i++)
160  if (partitioning[i] == MyRank)
161  {
162  mesh.GetElementVertices(i, vert);
163  element_counter++;
164  for (j = 0; j < vert.Size(); j++)
165  if (vert_global_local[vert[j]] < 0)
166  {
167  vert_global_local[vert[j]] = vert_counter++;
168  }
169  }
170 
171  NumOfVertices = vert_counter;
172  NumOfElements = element_counter;
173  vertices.SetSize(NumOfVertices);
174 
175  // re-enumerate the local vertices to preserve the global ordering
176  for (i = vert_counter = 0; i < vert_global_local.Size(); i++)
177  if (vert_global_local[i] >= 0)
178  {
179  vert_global_local[i] = vert_counter++;
180  }
181 
182  // determine vertices
183  for (i = 0; i < vert_global_local.Size(); i++)
184  if (vert_global_local[i] >= 0)
185  {
186  vertices[vert_global_local[i]].SetCoords(mesh.GetVertex(i));
187  }
188 
189  // determine elements
190  element_counter = 0;
191  elements.SetSize(NumOfElements);
192  for (i = 0; i < mesh.GetNE(); i++)
193  if (partitioning[i] == MyRank)
194  {
195  elements[element_counter] = mesh.GetElement(i)->Duplicate(this);
196  int *v = elements[element_counter]->GetVertices();
197  int nv = elements[element_counter]->GetNVertices();
198  for (j = 0; j < nv; j++)
199  {
200  v[j] = vert_global_local[v[j]];
201  }
202  element_counter++;
203  }
204 
205  Table *edge_element = NULL;
206  if (mesh.NURBSext)
207  {
208  activeBdrElem.SetSize(mesh.GetNBE());
209  activeBdrElem = false;
210  }
211  // build boundary elements
212  if (Dim == 3)
213  {
214  NumOfBdrElements = 0;
215  for (i = 0; i < mesh.GetNBE(); i++)
216  {
217  int face, o, el1, el2;
218  mesh.GetBdrElementFace(i, &face, &o);
219  mesh.GetFaceElements(face, &el1, &el2);
220  if (partitioning[(o % 2 == 0 || el2 < 0) ? el1 : el2] == MyRank)
221  {
223  if (mesh.NURBSext)
224  {
225  activeBdrElem[i] = true;
226  }
227  }
228  }
229 
230  bdrelem_counter = 0;
231  boundary.SetSize(NumOfBdrElements);
232  for (i = 0; i < mesh.GetNBE(); i++)
233  {
234  int face, o, el1, el2;
235  mesh.GetBdrElementFace(i, &face, &o);
236  mesh.GetFaceElements(face, &el1, &el2);
237  if (partitioning[(o % 2 == 0 || el2 < 0) ? el1 : el2] == MyRank)
238  {
239  boundary[bdrelem_counter] = mesh.GetBdrElement(i)->Duplicate(this);
240  int *v = boundary[bdrelem_counter]->GetVertices();
241  int nv = boundary[bdrelem_counter]->GetNVertices();
242  for (j = 0; j < nv; j++)
243  {
244  v[j] = vert_global_local[v[j]];
245  }
246  bdrelem_counter++;
247  }
248  }
249  }
250  else if (Dim == 2)
251  {
252  edge_element = new Table;
253  Transpose(mesh.ElementToEdgeTable(), *edge_element, mesh.GetNEdges());
254 
255  NumOfBdrElements = 0;
256  for (i = 0; i < mesh.GetNBE(); i++)
257  {
258  int edge = mesh.GetBdrElementEdgeIndex(i);
259  int el1 = edge_element->GetRow(edge)[0];
260  if (partitioning[el1] == MyRank)
261  {
262  NumOfBdrElements++;
263  if (mesh.NURBSext)
264  {
265  activeBdrElem[i] = true;
266  }
267  }
268  }
269 
270  bdrelem_counter = 0;
271  boundary.SetSize(NumOfBdrElements);
272  for (i = 0; i < mesh.GetNBE(); i++)
273  {
274  int edge = mesh.GetBdrElementEdgeIndex(i);
275  int el1 = edge_element->GetRow(edge)[0];
276  if (partitioning[el1] == MyRank)
277  {
278  boundary[bdrelem_counter] = mesh.GetBdrElement(i)->Duplicate(this);
279  int *v = boundary[bdrelem_counter]->GetVertices();
280  int nv = boundary[bdrelem_counter]->GetNVertices();
281  for (j = 0; j < nv; j++)
282  {
283  v[j] = vert_global_local[v[j]];
284  }
285  bdrelem_counter++;
286  }
287  }
288  }
289  else if (Dim == 1)
290  {
291  NumOfBdrElements = 0;
292  for (i = 0; i < mesh.GetNBE(); i++)
293  {
294  int vert = mesh.boundary[i]->GetVertices()[0];
295  int el1, el2;
296  mesh.GetFaceElements(vert, &el1, &el2);
297  if (partitioning[el1] == MyRank)
298  {
300  }
301  }
302 
303  bdrelem_counter = 0;
304  boundary.SetSize(NumOfBdrElements);
305  for (i = 0; i < mesh.GetNBE(); i++)
306  {
307  int vert = mesh.boundary[i]->GetVertices()[0];
308  int el1, el2;
309  mesh.GetFaceElements(vert, &el1, &el2);
310  if (partitioning[el1] == MyRank)
311  {
312  boundary[bdrelem_counter] = mesh.GetBdrElement(i)->Duplicate(this);
313  int *v = boundary[bdrelem_counter]->GetVertices();
314  v[0] = vert_global_local[v[0]];
315  bdrelem_counter++;
316  }
317  }
318  }
319 
320  meshgen = mesh.MeshGenerator();
321 
322  mesh.attributes.Copy(attributes);
324 
325  // this is called by the default Mesh constructor
326  // InitTables();
327 
328  if (Dim > 1)
329  {
330  el_to_edge = new Table;
332  }
333  else
334  {
335  NumOfEdges = 0;
336  }
337 
338  STable3D *faces_tbl = NULL;
339  if (Dim == 3)
340  {
341  faces_tbl = GetElementToFaceTable(1);
342  }
343  else
344  {
345  NumOfFaces = 0;
346  }
347  GenerateFaces();
348 
349  ListOfIntegerSets groups;
350  IntegerSet group;
351 
352  // the first group is the local one
353  group.Recreate(1, &MyRank);
354  groups.Insert(group);
355 
356 #ifdef MFEM_DEBUG
357  if (Dim < 3 && mesh.GetNFaces() != 0)
358  {
359  cerr << "ParMesh::ParMesh (proc " << MyRank << ") : "
360  "(Dim < 3 && mesh.GetNFaces() != 0) is true!" << endl;
361  mfem_error();
362  }
363 #endif
364  // determine shared faces
365  int sface_counter = 0;
366  Array<int> face_group(mesh.GetNFaces());
367  for (i = 0; i < face_group.Size(); i++)
368  {
369  int el[2];
370  face_group[i] = -1;
371  mesh.GetFaceElements(i, &el[0], &el[1]);
372  if (el[1] >= 0)
373  {
374  el[0] = partitioning[el[0]];
375  el[1] = partitioning[el[1]];
376  if ((el[0] == MyRank && el[1] != MyRank) ||
377  (el[0] != MyRank && el[1] == MyRank))
378  {
379  group.Recreate(2, el);
380  face_group[i] = groups.Insert(group) - 1;
381  sface_counter++;
382  }
383  }
384  }
385 
386  // determine shared edges
387  int sedge_counter = 0;
388  if (!edge_element)
389  {
390  edge_element = new Table;
391  if (Dim == 1)
392  {
393  edge_element->SetDims(0,0);
394  }
395  else
396  {
397  Transpose(mesh.ElementToEdgeTable(), *edge_element, mesh.GetNEdges());
398  }
399  }
400  for (i = 0; i < edge_element->Size(); i++)
401  {
402  int me = 0, others = 0;
403  for (j = edge_element->GetI()[i]; j < edge_element->GetI()[i+1]; j++)
404  {
405  edge_element->GetJ()[j] = partitioning[edge_element->GetJ()[j]];
406  if (edge_element->GetJ()[j] == MyRank)
407  {
408  me = 1;
409  }
410  else
411  {
412  others = 1;
413  }
414  }
415 
416  if (me && others)
417  {
418  sedge_counter++;
419  group.Recreate(edge_element->RowSize(i), edge_element->GetRow(i));
420  edge_element->GetRow(i)[0] = groups.Insert(group) - 1;
421  }
422  else
423  {
424  edge_element->GetRow(i)[0] = -1;
425  }
426  }
427 
428  // determine shared vertices
429  int svert_counter = 0;
430  Table *vert_element = mesh.GetVertexToElementTable(); // we must delete this
431 
432  for (i = 0; i < vert_element->Size(); i++)
433  {
434  int me = 0, others = 0;
435  for (j = vert_element->GetI()[i]; j < vert_element->GetI()[i+1]; j++)
436  {
437  vert_element->GetJ()[j] = partitioning[vert_element->GetJ()[j]];
438  if (vert_element->GetJ()[j] == MyRank)
439  {
440  me = 1;
441  }
442  else
443  {
444  others = 1;
445  }
446  }
447 
448  if (me && others)
449  {
450  svert_counter++;
451  group.Recreate(vert_element->RowSize(i), vert_element->GetRow(i));
452  vert_element->GetI()[i] = groups.Insert(group) - 1;
453  }
454  else
455  {
456  vert_element->GetI()[i] = -1;
457  }
458  }
459 
460  // build group_sface
461  group_sface.MakeI(groups.Size()-1);
462 
463  for (i = 0; i < face_group.Size(); i++)
464  if (face_group[i] >= 0)
465  {
466  group_sface.AddAColumnInRow(face_group[i]);
467  }
468 
469  group_sface.MakeJ();
470 
471  sface_counter = 0;
472  for (i = 0; i < face_group.Size(); i++)
473  if (face_group[i] >= 0)
474  {
475  group_sface.AddConnection(face_group[i], sface_counter++);
476  }
477 
479 
480  // build group_sedge
481  group_sedge.MakeI(groups.Size()-1);
482 
483  for (i = 0; i < edge_element->Size(); i++)
484  if (edge_element->GetRow(i)[0] >= 0)
485  {
486  group_sedge.AddAColumnInRow(edge_element->GetRow(i)[0]);
487  }
488 
489  group_sedge.MakeJ();
490 
491  sedge_counter = 0;
492  for (i = 0; i < edge_element->Size(); i++)
493  if (edge_element->GetRow(i)[0] >= 0)
494  group_sedge.AddConnection(edge_element->GetRow(i)[0],
495  sedge_counter++);
496 
498 
499  // build group_svert
500  group_svert.MakeI(groups.Size()-1);
501 
502  for (i = 0; i < vert_element->Size(); i++)
503  if (vert_element->GetI()[i] >= 0)
504  {
505  group_svert.AddAColumnInRow(vert_element->GetI()[i]);
506  }
507 
508  group_svert.MakeJ();
509 
510  svert_counter = 0;
511  for (i = 0; i < vert_element->Size(); i++)
512  if (vert_element->GetI()[i] >= 0)
513  group_svert.AddConnection(vert_element->GetI()[i],
514  svert_counter++);
515 
517 
518  // build shared_faces and sface_lface
519  shared_faces.SetSize(sface_counter);
520  sface_lface. SetSize(sface_counter);
521 
522  if (Dim == 3)
523  {
524  sface_counter = 0;
525  for (i = 0; i < face_group.Size(); i++)
526  if (face_group[i] >= 0)
527  {
528  shared_faces[sface_counter] = mesh.GetFace(i)->Duplicate(this);
529  int *v = shared_faces[sface_counter]->GetVertices();
530  int nv = shared_faces[sface_counter]->GetNVertices();
531  for (j = 0; j < nv; j++)
532  {
533  v[j] = vert_global_local[v[j]];
534  }
535  switch (shared_faces[sface_counter]->GetType())
536  {
537  case Element::TRIANGLE:
538  sface_lface[sface_counter] = (*faces_tbl)(v[0], v[1], v[2]);
539  // mark the shared face for refinement by reorienting
540  // it according to the refinement flag in the tetrahedron
541  // to which this shared face belongs to.
542  {
543  int lface = sface_lface[sface_counter];
544  Tetrahedron *tet =
545  (Tetrahedron *)(elements[faces_info[lface].Elem1No]);
546  int re[2], type, flag, *tv;
547  tet->ParseRefinementFlag(re, type, flag);
548  tv = tet->GetVertices();
549  switch (faces_info[lface].Elem1Inf/64)
550  {
551  case 0:
552  switch (re[1])
553  {
554  case 1: v[0] = tv[1]; v[1] = tv[2]; v[2] = tv[3];
555  break;
556  case 4: v[0] = tv[3]; v[1] = tv[1]; v[2] = tv[2];
557  break;
558  case 5: v[0] = tv[2]; v[1] = tv[3]; v[2] = tv[1];
559  break;
560  }
561  break;
562  case 1:
563  switch (re[0])
564  {
565  case 2: v[0] = tv[2]; v[1] = tv[0]; v[2] = tv[3];
566  break;
567  case 3: v[0] = tv[0]; v[1] = tv[3]; v[2] = tv[2];
568  break;
569  case 5: v[0] = tv[3]; v[1] = tv[2]; v[2] = tv[0];
570  break;
571  }
572  break;
573  case 2:
574  v[0] = tv[0]; v[1] = tv[1]; v[2] = tv[3];
575  break;
576  case 3:
577  v[0] = tv[1]; v[1] = tv[0]; v[2] = tv[2];
578  break;
579  }
580  // flip the shared face in the processor that owns the
581  // second element (in 'mesh')
582  {
583  int gl_el1, gl_el2;
584  mesh.GetFaceElements(i, &gl_el1, &gl_el2);
585  if (MyRank == partitioning[gl_el2])
586  {
587  const int t = v[0]; v[0] = v[1]; v[1] = t;
588  }
589  }
590  }
591  break;
593  sface_lface[sface_counter] =
594  (*faces_tbl)(v[0], v[1], v[2], v[3]);
595  break;
596  }
597  sface_counter++;
598  }
599 
600  delete faces_tbl;
601  }
602 
603  // build shared_edges and sedge_ledge
604  shared_edges.SetSize(sedge_counter);
605  sedge_ledge. SetSize(sedge_counter);
606 
607  {
608  DSTable v_to_v(NumOfVertices);
609  GetVertexToVertexTable(v_to_v);
610 
611  sedge_counter = 0;
612  for (i = 0; i < edge_element->Size(); i++)
613  if (edge_element->GetRow(i)[0] >= 0)
614  {
615  mesh.GetEdgeVertices(i, vert);
616 
617  shared_edges[sedge_counter] =
618  new Segment(vert_global_local[vert[0]],
619  vert_global_local[vert[1]], 1);
620 
621  if ((sedge_ledge[sedge_counter] =
622  v_to_v(vert_global_local[vert[0]],
623  vert_global_local[vert[1]])) < 0)
624  {
625  cerr << "\n\n\n" << MyRank << ": ParMesh::ParMesh: "
626  << "ERROR in v_to_v\n\n" << endl;
627  mfem_error();
628  }
629 
630  sedge_counter++;
631  }
632  }
633 
634  delete edge_element;
635 
636  // build svert_lvert
637  svert_lvert.SetSize(svert_counter);
638 
639  svert_counter = 0;
640  for (i = 0; i < vert_element->Size(); i++)
641  if (vert_element->GetI()[i] >= 0)
642  {
643  svert_lvert[svert_counter++] = vert_global_local[i];
644  }
645 
646  delete vert_element;
647 
648  // build the group communication topology
649  gtopo.Create(groups, 822);
650 
651  if (mesh.NURBSext)
652  {
653  NURBSext = new ParNURBSExtension(comm, mesh.NURBSext, partitioning,
654  activeBdrElem);
655  }
656 
657  if (mesh.GetNodes()) // curved mesh
658  {
659  Nodes = new ParGridFunction(this, mesh.GetNodes());
660  own_nodes = 1;
661 
662  Array<int> gvdofs, lvdofs;
663  Vector lnodes;
664  element_counter = 0;
665  for (i = 0; i < mesh.GetNE(); i++)
666  if (partitioning[i] == MyRank)
667  {
668  Nodes->FESpace()->GetElementVDofs(element_counter, lvdofs);
669  mesh.GetNodes()->FESpace()->GetElementVDofs(i, gvdofs);
670  mesh.GetNodes()->GetSubVector(gvdofs, lnodes);
671  Nodes->SetSubVector(lvdofs, lnodes);
672  element_counter++;
673  }
674  }
675 
676  if (partitioning_ == NULL)
677  {
678  delete [] partitioning;
679  }
680 
681  have_face_nbr_data = false;
682 }
683 
685  : MyComm(pncmesh.MyComm)
686  , NRanks(pncmesh.NRanks)
687  , MyRank(pncmesh.MyRank)
688  , gtopo(MyComm)
689  , pncmesh(NULL)
690 {
691  Mesh::Init();
693  Mesh::InitFromNCMesh(pncmesh);
694  have_face_nbr_data = false;
695 }
696 
697 void ParMesh::GroupEdge(int group, int i, int &edge, int &o)
698 {
699  int sedge = group_sedge.GetJ()[group_sedge.GetI()[group-1]+i];
700  edge = sedge_ledge[sedge];
701  int *v = shared_edges[sedge]->GetVertices();
702  o = (v[0] < v[1]) ? (+1) : (-1);
703 }
704 
705 void ParMesh::GroupFace(int group, int i, int &face, int &o)
706 {
707  int sface = group_sface.GetJ()[group_sface.GetI()[group-1]+i];
708  face = sface_lface[sface];
709  // face gives the base orientation
710  if (faces[face]->GetType() == Element::TRIANGLE)
711  o = GetTriOrientation(faces[face]->GetVertices(),
712  shared_faces[sface]->GetVertices());
713  if (faces[face]->GetType() == Element::QUADRILATERAL)
714  o = GetQuadOrientation(faces[face]->GetVertices(),
715  shared_faces[sface]->GetVertices());
716 }
717 
718 // For a line segment with vertices v[0] and v[1], return a number with
719 // the following meaning:
720 // 0 - the edge was not refined
721 // 1 - the edge e was refined once by splitting v[0],v[1]
722 int ParMesh::GetEdgeSplittings(Element *edge, const DSTable &v_to_v,
723  int *middle)
724 {
725  int m, *v = edge->GetVertices();
726 
727  if ((m = v_to_v(v[0], v[1])) != -1 && middle[m] != -1)
728  {
729  return 1;
730  }
731  else
732  {
733  return 0;
734  }
735 }
736 
737 // For a triangular face with (correctly ordered) vertices v[0], v[1], v[2]
738 // return a number with the following meaning:
739 // 0 - the face was not refined
740 // 1 - the face was refined once by splitting v[0],v[1]
741 // 2 - the face was refined twice by splitting v[0],v[1] and then v[1],v[2]
742 // 3 - the face was refined twice by splitting v[0],v[1] and then v[0],v[2]
743 // 4 - the face was refined three times (as in 2+3)
744 int ParMesh::GetFaceSplittings(Element *face, const DSTable &v_to_v,
745  int *middle)
746 {
747  int m, right = 0;
748  int number_of_splittings = 0;
749  int *v = face->GetVertices();
750 
751  if ((m = v_to_v(v[0], v[1])) != -1 && middle[m] != -1)
752  {
753  number_of_splittings++;
754  if ((m = v_to_v(v[1], v[2])) != -1 && middle[m] != -1)
755  {
756  right = 1;
757  number_of_splittings++;
758  }
759  if ((m = v_to_v(v[2], v[0])) != -1 && middle[m] != -1)
760  {
761  number_of_splittings++;
762  }
763 
764  switch (number_of_splittings)
765  {
766  case 2:
767  if (right == 0)
768  {
769  number_of_splittings++;
770  }
771  break;
772  case 3:
773  number_of_splittings++;
774  break;
775  }
776  }
777 
778  return number_of_splittings;
779 }
780 
781 void ParMesh::GenerateOffsets(int N, HYPRE_Int loc_sizes[],
782  Array<HYPRE_Int> *offsets[]) const
783 {
784  if (HYPRE_AssumedPartitionCheck())
785  {
786  Array<HYPRE_Int> temp(N);
787  MPI_Scan(loc_sizes, temp.GetData(), N, HYPRE_MPI_INT, MPI_SUM, MyComm);
788  for (int i = 0; i < N; i++)
789  {
790  offsets[i]->SetSize(3);
791  (*offsets[i])[0] = temp[i] - loc_sizes[i];
792  (*offsets[i])[1] = temp[i];
793  }
794  MPI_Bcast(temp.GetData(), N, HYPRE_MPI_INT, NRanks-1, MyComm);
795  for (int i = 0; i < N; i++)
796  {
797  (*offsets[i])[2] = temp[i];
798  // check for overflow
799  MFEM_VERIFY((*offsets[i])[0] >= 0 && (*offsets[i])[1] >= 0,
800  "overflow in offsets");
801  }
802  }
803  else
804  {
805  Array<HYPRE_Int> temp(N*NRanks);
806  MPI_Allgather(loc_sizes, N, HYPRE_MPI_INT, temp.GetData(), N,
807  HYPRE_MPI_INT, MyComm);
808  for (int i = 0; i < N; i++)
809  {
810  Array<HYPRE_Int> &offs = *offsets[i];
811  offs.SetSize(NRanks+1);
812  offs[0] = 0;
813  for (int j = 0; j < NRanks; j++)
814  {
815  offs[j+1] = offs[j] + temp[i+N*j];
816  }
817  // Check for overflow
818  MFEM_VERIFY(offs[MyRank] >= 0 && offs[MyRank+1] >= 0,
819  "overflow in offsets");
820  }
821  }
822 }
823 
825  int i, IsoparametricTransformation *ElTr)
826 {
827  DenseMatrix &pointmat = ElTr->GetPointMat();
828  Element *elem = face_nbr_elements[i];
829 
830  ElTr->Attribute = elem->GetAttribute();
831  ElTr->ElementNo = NumOfElements + i;
832 
833  if (Nodes == NULL)
834  {
835  const int nv = elem->GetNVertices();
836  const int *v = elem->GetVertices();
837 
838  pointmat.SetSize(spaceDim, nv);
839  for (int k = 0; k < spaceDim; k++)
840  {
841  for (int j = 0; j < nv; j++)
842  {
843  pointmat(k, j) = face_nbr_vertices[v[j]](k);
844  }
845  }
846 
848  }
849  else
850  {
851  Array<int> vdofs;
852  ParGridFunction *pNodes = dynamic_cast<ParGridFunction *>(Nodes);
853  if (pNodes)
854  {
855  pNodes->ParFESpace()->GetFaceNbrElementVDofs(i, vdofs);
856  int n = vdofs.Size()/spaceDim;
857  pointmat.SetSize(spaceDim, n);
858  for (int k = 0; k < spaceDim; k++)
859  {
860  for (int j = 0; j < n; j++)
861  {
862  pointmat(k,j) = (pNodes->FaceNbrData())(vdofs[n*k+j]);
863  }
864  }
865 
866  ElTr->SetFE(pNodes->ParFESpace()->GetFaceNbrFE(i));
867  }
868  else
869  {
870  MFEM_ABORT("Nodes are not ParGridFunction!");
871  }
872  }
873 }
874 
876 {
877  if (!have_face_nbr_data)
878  {
879  return;
880  }
881 
882  have_face_nbr_data = false;
886  for (int i = 0; i < face_nbr_elements.Size(); i++)
887  {
889  }
890  face_nbr_elements.DeleteAll();
891  face_nbr_vertices.DeleteAll();
894 }
895 
897 {
898  if (have_face_nbr_data)
899  {
900  return;
901  }
902 
903  if (Nonconforming())
904  {
905  // with ParNCMesh we can set up face neighbors without communication
906  pncmesh->GetFaceNeighbors(*this);
907  have_face_nbr_data = true;
908 
910  return;
911  }
912 
913  Table *gr_sface;
914  int *s2l_face;
915  if (Dim == 1)
916  {
917  gr_sface = &group_svert;
918  s2l_face = svert_lvert;
919  }
920  else if (Dim == 2)
921  {
922  gr_sface = &group_sedge;
923  s2l_face = sedge_ledge;
924  }
925  else
926  {
927  gr_sface = &group_sface;
928  s2l_face = sface_lface;
929  }
930 
931  int num_face_nbrs = 0;
932  for (int g = 1; g < GetNGroups(); g++)
933  if (gr_sface->RowSize(g-1) > 0)
934  {
935  num_face_nbrs++;
936  }
937 
938  face_nbr_group.SetSize(num_face_nbrs);
939 
940  if (num_face_nbrs == 0)
941  {
942  have_face_nbr_data = true;
943  return;
944  }
945 
946  {
947  // sort face-neighbors by processor rank
948  Array<Pair<int, int> > rank_group(num_face_nbrs);
949 
950  for (int g = 1, counter = 0; g < GetNGroups(); g++)
951  if (gr_sface->RowSize(g-1) > 0)
952  {
953 #ifdef MFEM_DEBUG
954  if (gtopo.GetGroupSize(g) != 2)
955  mfem_error("ParMesh::ExchangeFaceNbrData() : "
956  "group size is not 2!");
957 #endif
958  const int *nbs = gtopo.GetGroup(g);
959  int lproc = (nbs[0]) ? nbs[0] : nbs[1];
960  rank_group[counter].one = gtopo.GetNeighborRank(lproc);
961  rank_group[counter].two = g;
962  counter++;
963  }
964 
965  SortPairs<int, int>(rank_group, rank_group.Size());
966 
967  for (int fn = 0; fn < num_face_nbrs; fn++)
968  {
969  face_nbr_group[fn] = rank_group[fn].two;
970  }
971  }
972 
973  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
974  MPI_Request *send_requests = requests;
975  MPI_Request *recv_requests = requests + num_face_nbrs;
976  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
977 
978  int *nbr_data = new int[6*num_face_nbrs];
979  int *nbr_send_data = nbr_data;
980  int *nbr_recv_data = nbr_data + 3*num_face_nbrs;
981 
982  Array<int> el_marker(GetNE());
983  Array<int> vertex_marker(GetNV());
984  el_marker = -1;
985  vertex_marker = -1;
986 
987  Table send_face_nbr_elemdata, send_face_nbr_facedata;
988 
989  send_face_nbr_elements.MakeI(num_face_nbrs);
990  send_face_nbr_vertices.MakeI(num_face_nbrs);
991  send_face_nbr_elemdata.MakeI(num_face_nbrs);
992  send_face_nbr_facedata.MakeI(num_face_nbrs);
993  for (int fn = 0; fn < num_face_nbrs; fn++)
994  {
995  int nbr_group = face_nbr_group[fn];
996  int num_sfaces = gr_sface->RowSize(nbr_group-1);
997  int *sface = gr_sface->GetRow(nbr_group-1);
998  for (int i = 0; i < num_sfaces; i++)
999  {
1000  int lface = s2l_face[sface[i]];
1001  int el = faces_info[lface].Elem1No;
1002  if (el_marker[el] != fn)
1003  {
1004  el_marker[el] = fn;
1006 
1007  const int nv = elements[el]->GetNVertices();
1008  const int *v = elements[el]->GetVertices();
1009  for (int j = 0; j < nv; j++)
1010  if (vertex_marker[v[j]] != fn)
1011  {
1012  vertex_marker[v[j]] = fn;
1014  }
1015 
1016  send_face_nbr_elemdata.AddColumnsInRow(fn, nv + 2);
1017  }
1018  }
1019  send_face_nbr_facedata.AddColumnsInRow(fn, 2*num_sfaces);
1020 
1021  nbr_send_data[3*fn ] = send_face_nbr_elements.GetI()[fn];
1022  nbr_send_data[3*fn+1] = send_face_nbr_vertices.GetI()[fn];
1023  nbr_send_data[3*fn+2] = send_face_nbr_elemdata.GetI()[fn];
1024 
1025  int nbr_rank = GetFaceNbrRank(fn);
1026  int tag = 0;
1027 
1028  MPI_Isend(&nbr_send_data[3*fn], 3, MPI_INT, nbr_rank, tag, MyComm,
1029  &send_requests[fn]);
1030  MPI_Irecv(&nbr_recv_data[3*fn], 3, MPI_INT, nbr_rank, tag, MyComm,
1031  &recv_requests[fn]);
1032  }
1035  send_face_nbr_elemdata.MakeJ();
1036  send_face_nbr_facedata.MakeJ();
1037  el_marker = -1;
1038  vertex_marker = -1;
1039  for (int fn = 0; fn < num_face_nbrs; fn++)
1040  {
1041  int nbr_group = face_nbr_group[fn];
1042  int num_sfaces = gr_sface->RowSize(nbr_group-1);
1043  int *sface = gr_sface->GetRow(nbr_group-1);
1044  for (int i = 0; i < num_sfaces; i++)
1045  {
1046  int lface = s2l_face[sface[i]];
1047  int el = faces_info[lface].Elem1No;
1048  if (el_marker[el] != fn)
1049  {
1050  el_marker[el] = fn;
1052 
1053  const int nv = elements[el]->GetNVertices();
1054  const int *v = elements[el]->GetVertices();
1055  for (int j = 0; j < nv; j++)
1056  if (vertex_marker[v[j]] != fn)
1057  {
1058  vertex_marker[v[j]] = fn;
1060  }
1061 
1062  send_face_nbr_elemdata.AddConnection(fn, GetAttribute(el));
1063  send_face_nbr_elemdata.AddConnection(
1064  fn, GetElementBaseGeometry(el));
1065  send_face_nbr_elemdata.AddConnections(fn, v, nv);
1066  }
1067  send_face_nbr_facedata.AddConnection(fn, el);
1068  int info = faces_info[lface].Elem1Inf;
1069  // change the orientation in info to be relative to the shared face
1070  // in 1D and 2D keep the orientation equal to 0
1071  if (Dim == 3)
1072  {
1073  Element *lf = faces[lface];
1074  const int *sf_v = shared_faces[sface[i]]->GetVertices();
1075 
1076  if (lf->GetGeometryType() == Geometry::TRIANGLE)
1077  {
1078  info += GetTriOrientation(sf_v, lf->GetVertices());
1079  }
1080  else
1081  {
1082  info += GetQuadOrientation(sf_v, lf->GetVertices());
1083  }
1084  }
1085  send_face_nbr_facedata.AddConnection(fn, info);
1086  }
1087  }
1090  send_face_nbr_elemdata.ShiftUpI();
1091  send_face_nbr_facedata.ShiftUpI();
1092 
1093  // convert the vertex indices in send_face_nbr_elemdata
1094  // convert the element indices in send_face_nbr_facedata
1095  for (int fn = 0; fn < num_face_nbrs; fn++)
1096  {
1097  int num_elems = send_face_nbr_elements.RowSize(fn);
1098  int *elems = send_face_nbr_elements.GetRow(fn);
1099  int num_verts = send_face_nbr_vertices.RowSize(fn);
1100  int *verts = send_face_nbr_vertices.GetRow(fn);
1101  int *elemdata = send_face_nbr_elemdata.GetRow(fn);
1102  int num_sfaces = send_face_nbr_facedata.RowSize(fn)/2;
1103  int *facedata = send_face_nbr_facedata.GetRow(fn);
1104 
1105  for (int i = 0; i < num_verts; i++)
1106  {
1107  vertex_marker[verts[i]] = i;
1108  }
1109 
1110  for (int el = 0; el < num_elems; el++)
1111  {
1112  const int nv = elements[el]->GetNVertices();
1113  elemdata += 2; // skip the attribute and the geometry type
1114  for (int j = 0; j < nv; j++)
1115  {
1116  elemdata[j] = vertex_marker[elemdata[j]];
1117  }
1118  elemdata += nv;
1119 
1120  el_marker[elems[el]] = el;
1121  }
1122 
1123  for (int i = 0; i < num_sfaces; i++)
1124  {
1125  facedata[2*i] = el_marker[facedata[2*i]];
1126  }
1127  }
1128 
1129  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1130 
1131  Array<int> recv_face_nbr_facedata;
1132  Table recv_face_nbr_elemdata;
1133 
1134  // fill-in face_nbr_elements_offset, face_nbr_vertices_offset
1135  face_nbr_elements_offset.SetSize(num_face_nbrs + 1);
1136  face_nbr_vertices_offset.SetSize(num_face_nbrs + 1);
1137  recv_face_nbr_elemdata.MakeI(num_face_nbrs);
1138  face_nbr_elements_offset[0] = 0;
1139  face_nbr_vertices_offset[0] = 0;
1140  for (int fn = 0; fn < num_face_nbrs; fn++)
1141  {
1142  face_nbr_elements_offset[fn+1] =
1143  face_nbr_elements_offset[fn] + nbr_recv_data[3*fn];
1144  face_nbr_vertices_offset[fn+1] =
1145  face_nbr_vertices_offset[fn] + nbr_recv_data[3*fn+1];
1146  recv_face_nbr_elemdata.AddColumnsInRow(fn, nbr_recv_data[3*fn+2]);
1147  }
1148  recv_face_nbr_elemdata.MakeJ();
1149 
1150  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1151 
1152  // send and receive the element data
1153  for (int fn = 0; fn < num_face_nbrs; fn++)
1154  {
1155  int nbr_rank = GetFaceNbrRank(fn);
1156  int tag = 0;
1157 
1158  MPI_Isend(send_face_nbr_elemdata.GetRow(fn),
1159  send_face_nbr_elemdata.RowSize(fn),
1160  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
1161 
1162  MPI_Irecv(recv_face_nbr_elemdata.GetRow(fn),
1163  recv_face_nbr_elemdata.RowSize(fn),
1164  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
1165  }
1166 
1167  // convert the element data into face_nbr_elements
1168  face_nbr_elements.SetSize(face_nbr_elements_offset[num_face_nbrs]);
1169  while (true)
1170  {
1171  int fn;
1172  MPI_Waitany(num_face_nbrs, recv_requests, &fn, statuses);
1173 
1174  if (fn == MPI_UNDEFINED)
1175  {
1176  break;
1177  }
1178 
1179  int vert_off = face_nbr_vertices_offset[fn];
1180  int elem_off = face_nbr_elements_offset[fn];
1181  int num_elems = face_nbr_elements_offset[fn+1] - elem_off;
1182  int *recv_elemdata = recv_face_nbr_elemdata.GetRow(fn);
1183 
1184  for (int i = 0; i < num_elems; i++)
1185  {
1186  Element *el = NewElement(recv_elemdata[1]);
1187  el->SetAttribute(recv_elemdata[0]);
1188  recv_elemdata += 2;
1189  int nv = el->GetNVertices();
1190  for (int j = 0; j < nv; j++)
1191  {
1192  recv_elemdata[j] += vert_off;
1193  }
1194  el->SetVertices(recv_elemdata);
1195  recv_elemdata += nv;
1196  face_nbr_elements[elem_off++] = el;
1197  }
1198  }
1199 
1200  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1201 
1202  // send and receive the face data
1203  recv_face_nbr_facedata.SetSize(
1204  send_face_nbr_facedata.Size_of_connections());
1205  for (int fn = 0; fn < num_face_nbrs; fn++)
1206  {
1207  int nbr_rank = GetFaceNbrRank(fn);
1208  int tag = 0;
1209 
1210  MPI_Isend(send_face_nbr_facedata.GetRow(fn),
1211  send_face_nbr_facedata.RowSize(fn),
1212  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
1213 
1214  // the size of the send and receive face data is the same
1215  MPI_Irecv(&recv_face_nbr_facedata[send_face_nbr_facedata.GetI()[fn]],
1216  send_face_nbr_facedata.RowSize(fn),
1217  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
1218  }
1219 
1220  // transfer the received face data into faces_info
1221  while (true)
1222  {
1223  int fn;
1224  MPI_Waitany(num_face_nbrs, recv_requests, &fn, statuses);
1225 
1226  if (fn == MPI_UNDEFINED)
1227  {
1228  break;
1229  }
1230 
1231  int elem_off = face_nbr_elements_offset[fn];
1232  int nbr_group = face_nbr_group[fn];
1233  int num_sfaces = gr_sface->RowSize(nbr_group-1);
1234  int *sface = gr_sface->GetRow(nbr_group-1);
1235  int *facedata =
1236  &recv_face_nbr_facedata[send_face_nbr_facedata.GetI()[fn]];
1237 
1238  for (int i = 0; i < num_sfaces; i++)
1239  {
1240  int lface = s2l_face[sface[i]];
1241  FaceInfo &face_info = faces_info[lface];
1242  face_info.Elem2No = -1 - (facedata[2*i] + elem_off);
1243  int info = facedata[2*i+1];
1244  // change the orientation in info to be relative to the local face
1245  if (Dim < 3)
1246  {
1247  info++; // orientation 0 --> orientation 1
1248  }
1249  else
1250  {
1251  int nbr_ori = info%64, nbr_v[4];
1252  Element *lf = faces[lface];
1253  const int *sf_v = shared_faces[sface[i]]->GetVertices();
1254 
1255  if (lf->GetGeometryType() == Geometry::TRIANGLE)
1256  {
1257  // apply the nbr_ori to sf_v to get nbr_v
1258  const int *perm = tri_t::Orient[nbr_ori];
1259  for (int j = 0; j < 3; j++)
1260  {
1261  nbr_v[perm[j]] = sf_v[j];
1262  }
1263  // get the orientation of nbr_v w.r.t. the local face
1264  nbr_ori = GetTriOrientation(lf->GetVertices(), nbr_v);
1265  }
1266  else
1267  {
1268  // apply the nbr_ori to sf_v to get nbr_v
1269  const int *perm = quad_t::Orient[nbr_ori];
1270  for (int j = 0; j < 4; j++)
1271  {
1272  nbr_v[perm[j]] = sf_v[j];
1273  }
1274  // get the orientation of nbr_v w.r.t. the local face
1275  nbr_ori = GetQuadOrientation(lf->GetVertices(), nbr_v);
1276  }
1277 
1278  info = 64*(info/64) + nbr_ori;
1279  }
1280  face_info.Elem2Inf = info;
1281  }
1282  }
1283 
1284  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1285 
1286  // allocate the face_nbr_vertices
1287  face_nbr_vertices.SetSize(face_nbr_vertices_offset[num_face_nbrs]);
1288 
1289  delete [] nbr_data;
1290 
1291  delete [] statuses;
1292  delete [] requests;
1293 
1294  have_face_nbr_data = true;
1295 
1297 }
1298 
1300 {
1301  if (!have_face_nbr_data)
1302  {
1303  ExchangeFaceNbrData(); // calls this method at the end
1304  }
1305  else if (Nodes == NULL)
1306  {
1307  if (Nonconforming())
1308  {
1309  // with ParNCMesh we already have the vertices
1310  return;
1311  }
1312 
1313  int num_face_nbrs = GetNFaceNeighbors();
1314 
1315  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
1316  MPI_Request *send_requests = requests;
1317  MPI_Request *recv_requests = requests + num_face_nbrs;
1318  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
1319 
1320  // allocate buffer and copy the vertices to be sent
1322  for (int i = 0; i < send_vertices.Size(); i++)
1323  {
1324  send_vertices[i] = vertices[send_face_nbr_vertices.GetJ()[i]];
1325  }
1326 
1327  // send and receive the vertices
1328  for (int fn = 0; fn < num_face_nbrs; fn++)
1329  {
1330  int nbr_rank = GetFaceNbrRank(fn);
1331  int tag = 0;
1332 
1333  MPI_Isend(send_vertices[send_face_nbr_vertices.GetI()[fn]](),
1335  MPI_DOUBLE, nbr_rank, tag, MyComm, &send_requests[fn]);
1336 
1338  3*(face_nbr_vertices_offset[fn+1] -
1340  MPI_DOUBLE, nbr_rank, tag, MyComm, &recv_requests[fn]);
1341  }
1342 
1343  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1344  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1345 
1346  delete [] statuses;
1347  delete [] requests;
1348  }
1349  else
1350  {
1351  ParGridFunction *pNodes = dynamic_cast<ParGridFunction *>(Nodes);
1352  MFEM_VERIFY(pNodes != NULL, "Nodes are not ParGridFunction!");
1353  pNodes->ExchangeFaceNbrData();
1354  }
1355 }
1356 
1357 int ParMesh::GetFaceNbrRank(int fn) const
1358 {
1359  if (Conforming())
1360  {
1361  int nbr_group = face_nbr_group[fn];
1362  const int *nbs = gtopo.GetGroup(nbr_group);
1363  int nbr_lproc = (nbs[0]) ? nbs[0] : nbs[1];
1364  int nbr_rank = gtopo.GetNeighborRank(nbr_lproc);
1365  return nbr_rank;
1366  }
1367  else
1368  {
1369  // NC: simplified handling of face neighbor ranks
1370  return face_nbr_group[fn];
1371  }
1372 }
1373 
1375 {
1376  const Array<int> *s2l_face;
1377  if (Dim == 1)
1378  {
1379  s2l_face = &svert_lvert;
1380  }
1381  else if (Dim == 2)
1382  {
1383  s2l_face = &sedge_ledge;
1384  }
1385  else
1386  {
1387  s2l_face = &sface_lface;
1388  }
1389 
1390  Table *face_elem = new Table;
1391 
1392  face_elem->MakeI(faces_info.Size());
1393 
1394  for (int i = 0; i < faces_info.Size(); i++)
1395  {
1396  if (faces_info[i].Elem2No >= 0)
1397  {
1398  face_elem->AddColumnsInRow(i, 2);
1399  }
1400  else
1401  {
1402  face_elem->AddAColumnInRow(i);
1403  }
1404  }
1405  for (int i = 0; i < s2l_face->Size(); i++)
1406  {
1407  face_elem->AddAColumnInRow((*s2l_face)[i]);
1408  }
1409 
1410  face_elem->MakeJ();
1411 
1412  for (int i = 0; i < faces_info.Size(); i++)
1413  {
1414  face_elem->AddConnection(i, faces_info[i].Elem1No);
1415  if (faces_info[i].Elem2No >= 0)
1416  {
1417  face_elem->AddConnection(i, faces_info[i].Elem2No);
1418  }
1419  }
1420  for (int i = 0; i < s2l_face->Size(); i++)
1421  {
1422  int lface = (*s2l_face)[i];
1423  int nbr_elem_idx = -1 - faces_info[lface].Elem2No;
1424  face_elem->AddConnection(lface, NumOfElements + nbr_elem_idx);
1425  }
1426 
1427  face_elem->ShiftUpI();
1428 
1429  return face_elem;
1430 }
1431 
1433  FaceElementTransformations* FETr, int face_type, int face_geom)
1434 {
1435  // calculate composition of FETr->Loc1 and FETr->Elem1
1437  if (Nodes == NULL)
1438  {
1439  FETr->Elem1->Transform(FETr->Loc1.Transf.GetPointMat(), face_pm);
1441  }
1442  else
1443  {
1444  const FiniteElement* face_el =
1445  Nodes->FESpace()->GetTraceElement(FETr->Elem1No, face_geom);
1446 
1447 #if 0 // TODO: handle the case of non-interpolatory Nodes
1448  DenseMatrix I;
1449  face_el->Project(Transformation.GetFE(), FETr->Loc1.Transf, I);
1450  MultABt(Transformation.GetPointMat(), I, pm_face);
1451 #else
1452  IntegrationRule eir(face_el->GetDof());
1453  FETr->Loc1.Transform(face_el->GetNodes(), eir);
1454  Nodes->GetVectorValues(*FETr->Elem1, eir, face_pm);
1455 #endif
1456  FaceTransformation.SetFE(face_el);
1457  }
1458  return &FaceTransformation;
1459 }
1460 
1462 GetSharedFaceTransformations(int sf, bool fill2)
1463 {
1464  int FaceNo = GetSharedFace(sf);
1465 
1466  FaceInfo &face_info = faces_info[FaceNo];
1467 
1468  bool is_slave = Nonconforming() && IsSlaveFace(face_info);
1469  bool is_ghost = Nonconforming() && FaceNo >= GetNumFaces();
1470 
1471  NCFaceInfo* nc_info = NULL;
1472  if (is_slave) { nc_info = &nc_faces_info[face_info.NCFace]; }
1473 
1474  int local_face = is_ghost ? nc_info->MasterFace : FaceNo;
1475  int face_type = GetFaceElementType(local_face);
1476  int face_geom = GetFaceGeometryType(local_face);
1477 
1478  // setup the transformation for the first element
1479  FaceElemTr.Elem1No = face_info.Elem1No;
1482 
1483  // setup the transformation for the second (neighbor) element
1484  if (fill2)
1485  {
1486  FaceElemTr.Elem2No = -1 - face_info.Elem2No;
1489  }
1490  else
1491  {
1492  FaceElemTr.Elem2No = -1;
1493  }
1494 
1495  // setup the face transformation if the face is not a ghost
1496  FaceElemTr.FaceGeom = face_geom;
1497  if (!is_ghost)
1498  {
1500  // NOTE: The above call overwrites FaceElemTr.Loc1
1501  }
1502 
1503  // setup Loc1 & Loc2
1504  int elem_type = GetElementType(face_info.Elem1No);
1505  GetLocalFaceTransformation(face_type, elem_type, FaceElemTr.Loc1.Transf,
1506  face_info.Elem1Inf);
1507 
1508  if (fill2)
1509  {
1510  elem_type = face_nbr_elements[FaceElemTr.Elem2No]->GetType();
1511  GetLocalFaceTransformation(face_type, elem_type, FaceElemTr.Loc2.Transf,
1512  face_info.Elem2Inf);
1513  }
1514 
1515  // adjust Loc1 or Loc2 of the master face if this is a slave face
1516  if (is_slave)
1517  {
1518  // is a ghost slave? -> master not a ghost -> choose Elem1 local transf
1519  // not a ghost slave? -> master is a ghost -> choose Elem2 local transf
1521  is_ghost ? FaceElemTr.Loc1.Transf : FaceElemTr.Loc2.Transf;
1522 
1523  if (is_ghost || fill2)
1524  {
1525  ApplyLocalSlaveTransformation(loctr, face_info);
1526  }
1527 
1528  if (face_type == Element::SEGMENT && fill2)
1529  {
1530  // fix slave orientation in 2D: flip Loc2 to match Loc1 and Face
1532  std::swap(pm(0,0), pm(0,1));
1533  std::swap(pm(1,0), pm(1,1));
1534  }
1535  }
1536 
1537  // for ghost faces we need a special version of GetFaceTransformation
1538  if (is_ghost)
1539  {
1540  FaceElemTr.Face =
1541  GetGhostFaceTransformation(&FaceElemTr, face_type, face_geom);
1542  }
1543 
1544  return &FaceElemTr;
1545 }
1546 
1548 {
1549  if (Conforming())
1550  {
1551  switch (Dim)
1552  {
1553  case 1: return svert_lvert.Size();
1554  case 2: return sedge_ledge.Size();
1555  default: return sface_lface.Size();
1556  }
1557  }
1558  else
1559  {
1560  MFEM_ASSERT(Dim > 1, "");
1561  const NCMesh::NCList &shared = pncmesh->GetSharedList(Dim-1);
1562  return shared.conforming.size() + shared.slaves.size();
1563  }
1564 }
1565 
1566 int ParMesh::GetSharedFace(int sface) const
1567 {
1568  if (Conforming())
1569  {
1570  switch (Dim)
1571  {
1572  case 1: return svert_lvert[sface];
1573  case 2: return sedge_ledge[sface];
1574  default: return sface_lface[sface];
1575  }
1576  }
1577  else
1578  {
1579  MFEM_ASSERT(Dim > 1, "");
1580  const NCMesh::NCList &shared = pncmesh->GetSharedList(Dim-1);
1581  int csize = (int) shared.conforming.size();
1582  return sface < csize
1583  ? shared.conforming[sface].index
1584  : shared.slaves[sface - csize].index;
1585  }
1586 }
1587 
1589 {
1590  if (Dim != 3 || !(meshgen & 1))
1591  {
1592  return;
1593  }
1594 
1596 
1597  int *v;
1598 
1599  // The local edge and face numbering is changed therefore we need to
1600  // update sedge_ledge and sface_lface.
1601  {
1602  DSTable v_to_v(NumOfVertices);
1603  GetVertexToVertexTable(v_to_v);
1604  for (int i = 0; i < shared_edges.Size(); i++)
1605  {
1606  v = shared_edges[i]->GetVertices();
1607  sedge_ledge[i] = v_to_v(v[0], v[1]);
1608  }
1609  }
1610 
1611  // Rotate shared faces and update sface_lface.
1612  // Note that no communication is needed to ensure that the shared
1613  // faces are rotated in the same way in both processors. This is
1614  // automatic due to various things, e.g. the global to local vertex
1615  // mapping preserves the global order; also the way new vertices
1616  // are introduced during refinement is essential.
1617  {
1618  STable3D *faces_tbl = GetFacesTable();
1619  for (int i = 0; i < shared_faces.Size(); i++)
1620  if (shared_faces[i]->GetType() == Element::TRIANGLE)
1621  {
1622  v = shared_faces[i]->GetVertices();
1623 
1624  Rotate3(v[0], v[1], v[2]);
1625 
1626  sface_lface[i] = (*faces_tbl)(v[0], v[1], v[2]);
1627  }
1628  delete faces_tbl;
1629  }
1630 }
1631 
1632 void ParMesh::LocalRefinement(const Array<int> &marked_el, int type)
1633 {
1634  int i, j;
1635 
1636  if (pncmesh)
1637  {
1638  MFEM_ABORT("Local and nonconforming refinements cannot be mixed.");
1639  }
1640 
1642 
1644 
1645  if (Dim == 3)
1646  {
1647  int uniform_refinement = 0;
1648  if (type < 0)
1649  {
1650  type = -type;
1651  uniform_refinement = 1;
1652  }
1653 
1654  // 1. Get table of vertex to vertex connections.
1655  DSTable v_to_v(NumOfVertices);
1656  GetVertexToVertexTable(v_to_v);
1657 
1658  // 2. Create a marker array for all edges (vertex to vertex connections).
1659  Array<int> middle(v_to_v.NumberOfEntries());
1660  middle = -1;
1661 
1662  // 3. Do the red refinement.
1663  switch (type)
1664  {
1665  case 1:
1666  for (i = 0; i < marked_el.Size(); i++)
1667  {
1668  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1669  }
1670  break;
1671  case 2:
1672  for (i = 0; i < marked_el.Size(); i++)
1673  {
1674  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1675 
1676  Bisection(NumOfElements - 1, v_to_v, NULL, NULL, middle);
1677  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1678  }
1679  break;
1680  case 3:
1681  for (i = 0; i < marked_el.Size(); i++)
1682  {
1683  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1684 
1685  j = NumOfElements - 1;
1686  Bisection(j, v_to_v, NULL, NULL, middle);
1687  Bisection(NumOfElements - 1, v_to_v, NULL, NULL, middle);
1688  Bisection(j, v_to_v, NULL, NULL, middle);
1689 
1690  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1691  Bisection(NumOfElements-1, v_to_v, NULL, NULL, middle);
1692  Bisection(marked_el[i], v_to_v, NULL, NULL, middle);
1693  }
1694  break;
1695  }
1696 
1697  // 4. Do the green refinement (to get conforming mesh).
1698  int need_refinement;
1699  int refined_edge[5][3] =
1700  {
1701  {0, 0, 0},
1702  {1, 0, 0},
1703  {1, 1, 0},
1704  {1, 0, 1},
1705  {1, 1, 1}
1706  };
1707  int faces_in_group, max_faces_in_group = 0;
1708  // face_splittings identify how the shared faces have been split
1709  int **face_splittings = new int*[GetNGroups()-1];
1710  for (i = 0; i < GetNGroups()-1; i++)
1711  {
1712  faces_in_group = GroupNFaces(i+1);
1713  face_splittings[i] = new int[faces_in_group];
1714  if (faces_in_group > max_faces_in_group)
1715  {
1716  max_faces_in_group = faces_in_group;
1717  }
1718  }
1719  int neighbor, *iBuf = new int[max_faces_in_group];
1720 
1721  Array<int> group_faces;
1722 
1723  MPI_Request request;
1724  MPI_Status status;
1725 
1726 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1727  int ref_loops_all = 0, ref_loops_par = 0;
1728 #endif
1729  do
1730  {
1731  need_refinement = 0;
1732  for (i = 0; i < NumOfElements; i++)
1733  {
1734  if (elements[i]->NeedRefinement(v_to_v, middle))
1735  {
1736  need_refinement = 1;
1737  Bisection(i, v_to_v, NULL, NULL, middle);
1738  }
1739  }
1740 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1741  ref_loops_all++;
1742 #endif
1743 
1744  if (uniform_refinement)
1745  {
1746  continue;
1747  }
1748 
1749  // if the mesh is locally conforming start making it globally
1750  // conforming
1751  if (need_refinement == 0)
1752  {
1753 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1754  ref_loops_par++;
1755 #endif
1756  // MPI_Barrier(MyComm);
1757  const int tag = 293;
1758 
1759  // (a) send the type of interface splitting
1760  for (i = 0; i < GetNGroups()-1; i++)
1761  {
1762  group_sface.GetRow(i, group_faces);
1763  faces_in_group = group_faces.Size();
1764  // it is enough to communicate through the faces
1765  if (faces_in_group == 0) { continue; }
1766 
1767  for (j = 0; j < faces_in_group; j++)
1768  {
1769  face_splittings[i][j] =
1770  GetFaceSplittings(shared_faces[group_faces[j]], v_to_v,
1771  middle);
1772  }
1773  const int *nbs = gtopo.GetGroup(i+1);
1774  neighbor = gtopo.GetNeighborRank(nbs[0] ? nbs[0] : nbs[1]);
1775  MPI_Isend(face_splittings[i], faces_in_group, MPI_INT,
1776  neighbor, tag, MyComm, &request);
1777  }
1778 
1779  // (b) receive the type of interface splitting
1780  for (i = 0; i < GetNGroups()-1; i++)
1781  {
1782  group_sface.GetRow(i, group_faces);
1783  faces_in_group = group_faces.Size();
1784  if (faces_in_group == 0) { continue; }
1785 
1786  const int *nbs = gtopo.GetGroup(i+1);
1787  neighbor = gtopo.GetNeighborRank(nbs[0] ? nbs[0] : nbs[1]);
1788  MPI_Recv(iBuf, faces_in_group, MPI_INT, neighbor,
1789  tag, MyComm, &status);
1790 
1791  for (j = 0; j < faces_in_group; j++)
1792  {
1793  if (iBuf[j] == face_splittings[i][j]) { continue; }
1794 
1795  int *v = shared_faces[group_faces[j]]->GetVertices();
1796  for (int k = 0; k < 3; k++)
1797  {
1798  if (refined_edge[iBuf[j]][k] != 1 ||
1799  refined_edge[face_splittings[i][j]][k] != 0)
1800  { continue; }
1801 
1802  int ind[2] = { v[k], v[(k+1)%3] };
1803  int ii = v_to_v(ind[0], ind[1]);
1804  if (middle[ii] == -1)
1805  {
1806  need_refinement = 1;
1807  middle[ii] = NumOfVertices++;
1808  vertices.Append(Vertex());
1809  AverageVertices(ind, 2, vertices.Size()-1);
1810  }
1811  }
1812  }
1813  }
1814 
1815  i = need_refinement;
1816  MPI_Allreduce(&i, &need_refinement, 1, MPI_INT, MPI_LOR, MyComm);
1817  }
1818  }
1819  while (need_refinement == 1);
1820 
1821 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1822  i = ref_loops_all;
1823  MPI_Reduce(&i, &ref_loops_all, 1, MPI_INT, MPI_MAX, 0, MyComm);
1824  if (MyRank == 0)
1825  {
1826  cout << "\n\nParMesh::LocalRefinement : max. ref_loops_all = "
1827  << ref_loops_all << ", ref_loops_par = " << ref_loops_par
1828  << '\n' << endl;
1829  }
1830 #endif
1831 
1832  delete [] iBuf;
1833  for (i = 0; i < GetNGroups()-1; i++)
1834  {
1835  delete [] face_splittings[i];
1836  }
1837  delete [] face_splittings;
1838 
1839 
1840  // 5. Update the boundary elements.
1841  do
1842  {
1843  need_refinement = 0;
1844  for (i = 0; i < NumOfBdrElements; i++)
1845  {
1846  if (boundary[i]->NeedRefinement(v_to_v, middle))
1847  {
1848  need_refinement = 1;
1849  Bisection(i, v_to_v, middle);
1850  }
1851  }
1852  }
1853  while (need_refinement == 1);
1854 
1855  if (NumOfBdrElements != boundary.Size())
1856  mfem_error("ParMesh::LocalRefinement :"
1857  " (NumOfBdrElements != boundary.Size())");
1858 
1859  // 5a. Update the groups after refinement.
1860  if (el_to_face != NULL)
1861  {
1862  RefineGroups(v_to_v, middle);
1863  // GetElementToFaceTable(); // Called by RefineGroups
1864  GenerateFaces();
1865  }
1866 
1867  // 6. Un-mark the Pf elements.
1868  int refinement_edges[2], type, flag;
1869  for (i = 0; i < NumOfElements; i++)
1870  {
1871  Tetrahedron* el = (Tetrahedron*) elements[i];
1872  el->ParseRefinementFlag(refinement_edges, type, flag);
1873 
1874  if (type == Tetrahedron::TYPE_PF)
1875  {
1876  el->CreateRefinementFlag(refinement_edges, Tetrahedron::TYPE_PU,
1877  flag);
1878  }
1879  }
1880 
1881  // 7. Free the allocated memory.
1882  middle.DeleteAll();
1883 
1884  if (el_to_edge != NULL)
1885  {
1887  }
1888  } // 'if (Dim == 3)'
1889 
1890 
1891  if (Dim == 2)
1892  {
1893  int uniform_refinement = 0;
1894  if (type < 0)
1895  {
1896  type = -type;
1897  uniform_refinement = 1;
1898  }
1899 
1900  // 1. Get table of vertex to vertex connections.
1901  DSTable v_to_v(NumOfVertices);
1902  GetVertexToVertexTable(v_to_v);
1903 
1904  // 2. Get edge to element connections in arrays edge1 and edge2
1905  int nedges = v_to_v.NumberOfEntries();
1906  int *edge1 = new int[nedges];
1907  int *edge2 = new int[nedges];
1908  int *middle = new int[nedges];
1909 
1910  for (i = 0; i < nedges; i++)
1911  {
1912  edge1[i] = edge2[i] = middle[i] = -1;
1913  }
1914 
1915  for (i = 0; i < NumOfElements; i++)
1916  {
1917  int *v = elements[i]->GetVertices();
1918  for (j = 0; j < 3; j++)
1919  {
1920  int ind = v_to_v(v[j], v[(j+1)%3]);
1921  (edge1[ind] == -1) ? (edge1[ind] = i) : (edge2[ind] = i);
1922  }
1923  }
1924 
1925  // 3. Do the red refinement.
1926  for (i = 0; i < marked_el.Size(); i++)
1927  {
1928  RedRefinement(marked_el[i], v_to_v, edge1, edge2, middle);
1929  }
1930 
1931  // 4. Do the green refinement (to get conforming mesh).
1932  int need_refinement;
1933  int edges_in_group, max_edges_in_group = 0;
1934  // edge_splittings identify how the shared edges have been split
1935  int **edge_splittings = new int*[GetNGroups()-1];
1936  for (i = 0; i < GetNGroups()-1; i++)
1937  {
1938  edges_in_group = GroupNEdges(i+1);
1939  edge_splittings[i] = new int[edges_in_group];
1940  if (edges_in_group > max_edges_in_group)
1941  {
1942  max_edges_in_group = edges_in_group;
1943  }
1944  }
1945  int neighbor, *iBuf = new int[max_edges_in_group];
1946 
1947  Array<int> group_edges;
1948 
1949  MPI_Request request;
1950  MPI_Status status;
1951  Vertex V;
1952  V(2) = 0.0;
1953 
1954 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1955  int ref_loops_all = 0, ref_loops_par = 0;
1956 #endif
1957  do
1958  {
1959  need_refinement = 0;
1960  for (i = 0; i < nedges; i++)
1961  if (middle[i] != -1 && edge1[i] != -1)
1962  {
1963  need_refinement = 1;
1964  GreenRefinement(edge1[i], v_to_v, edge1, edge2, middle);
1965  }
1966 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1967  ref_loops_all++;
1968 #endif
1969 
1970  if (uniform_refinement)
1971  {
1972  continue;
1973  }
1974 
1975  // if the mesh is locally conforming start making it globally
1976  // conforming
1977  if (need_refinement == 0)
1978  {
1979 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
1980  ref_loops_par++;
1981 #endif
1982  // MPI_Barrier(MyComm);
1983 
1984  // (a) send the type of interface splitting
1985  for (i = 0; i < GetNGroups()-1; i++)
1986  {
1987  group_sedge.GetRow(i, group_edges);
1988  edges_in_group = group_edges.Size();
1989  // it is enough to communicate through the edges
1990  if (edges_in_group != 0)
1991  {
1992  for (j = 0; j < edges_in_group; j++)
1993  edge_splittings[i][j] =
1994  GetEdgeSplittings(shared_edges[group_edges[j]], v_to_v,
1995  middle);
1996  const int *nbs = gtopo.GetGroup(i+1);
1997  if (nbs[0] == 0)
1998  {
1999  neighbor = gtopo.GetNeighborRank(nbs[1]);
2000  }
2001  else
2002  {
2003  neighbor = gtopo.GetNeighborRank(nbs[0]);
2004  }
2005  MPI_Isend(edge_splittings[i], edges_in_group, MPI_INT,
2006  neighbor, 0, MyComm, &request);
2007  }
2008  }
2009 
2010  // (b) receive the type of interface splitting
2011  for (i = 0; i < GetNGroups()-1; i++)
2012  {
2013  group_sedge.GetRow(i, group_edges);
2014  edges_in_group = group_edges.Size();
2015  if (edges_in_group != 0)
2016  {
2017  const int *nbs = gtopo.GetGroup(i+1);
2018  if (nbs[0] == 0)
2019  {
2020  neighbor = gtopo.GetNeighborRank(nbs[1]);
2021  }
2022  else
2023  {
2024  neighbor = gtopo.GetNeighborRank(nbs[0]);
2025  }
2026  MPI_Recv(iBuf, edges_in_group, MPI_INT, neighbor,
2027  MPI_ANY_TAG, MyComm, &status);
2028 
2029  for (j = 0; j < edges_in_group; j++)
2030  if (iBuf[j] == 1 && edge_splittings[i][j] == 0)
2031  {
2032  int *v = shared_edges[group_edges[j]]->GetVertices();
2033  int ii = v_to_v(v[0], v[1]);
2034 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
2035  if (middle[ii] != -1)
2036  mfem_error("ParMesh::LocalRefinement (triangles) : "
2037  "Oops!");
2038 #endif
2039  need_refinement = 1;
2040  middle[ii] = NumOfVertices++;
2041  for (int c = 0; c < 2; c++)
2042  {
2043  V(c) = 0.5 * (vertices[v[0]](c) + vertices[v[1]](c));
2044  }
2045  vertices.Append(V);
2046  }
2047  }
2048  }
2049 
2050  i = need_refinement;
2051  MPI_Allreduce(&i, &need_refinement, 1, MPI_INT, MPI_LOR, MyComm);
2052  }
2053  }
2054  while (need_refinement == 1);
2055 
2056 #ifdef MFEM_DEBUG_PARMESH_LOCALREF
2057  i = ref_loops_all;
2058  MPI_Reduce(&i, &ref_loops_all, 1, MPI_INT, MPI_MAX, 0, MyComm);
2059  if (MyRank == 0)
2060  {
2061  cout << "\n\nParMesh::LocalRefinement : max. ref_loops_all = "
2062  << ref_loops_all << ", ref_loops_par = " << ref_loops_par
2063  << '\n' << endl;
2064  }
2065 #endif
2066 
2067  for (i = 0; i < GetNGroups()-1; i++)
2068  {
2069  delete [] edge_splittings[i];
2070  }
2071  delete [] edge_splittings;
2072 
2073  delete [] iBuf;
2074 
2075  // 5. Update the boundary elements.
2076  int v1[2], v2[2], bisect, temp;
2077  temp = NumOfBdrElements;
2078  for (i = 0; i < temp; i++)
2079  {
2080  int *v = boundary[i]->GetVertices();
2081  bisect = v_to_v(v[0], v[1]);
2082  if (middle[bisect] != -1)
2083  {
2084  // the element was refined (needs updating)
2085  if (boundary[i]->GetType() == Element::SEGMENT)
2086  {
2087  v1[0] = v[0]; v1[1] = middle[bisect];
2088  v2[0] = middle[bisect]; v2[1] = v[1];
2089 
2090  boundary[i]->SetVertices(v1);
2091  boundary.Append(new Segment(v2, boundary[i]->GetAttribute()));
2092  }
2093  else
2094  mfem_error("Only bisection of segment is implemented for bdr"
2095  " elem.");
2096  }
2097  }
2098  NumOfBdrElements = boundary.Size();
2099 
2100  // 5a. Update the groups after refinement.
2101  RefineGroups(v_to_v, middle);
2102 
2103  // 6. Free the allocated memory.
2104  delete [] edge1;
2105  delete [] edge2;
2106  delete [] middle;
2107 
2108  if (el_to_edge != NULL)
2109  {
2111  GenerateFaces();
2112  }
2113  } // 'if (Dim == 2)'
2114 
2115  if (Dim == 1) // --------------------------------------------------------
2116  {
2117  int cne = NumOfElements, cnv = NumOfVertices;
2118  NumOfVertices += marked_el.Size();
2119  NumOfElements += marked_el.Size();
2120  vertices.SetSize(NumOfVertices);
2121  elements.SetSize(NumOfElements);
2123 
2124  for (j = 0; j < marked_el.Size(); j++)
2125  {
2126  i = marked_el[j];
2127  Segment *c_seg = (Segment *)elements[i];
2128  int *vert = c_seg->GetVertices(), attr = c_seg->GetAttribute();
2129  int new_v = cnv + j, new_e = cne + j;
2130  AverageVertices(vert, 2, new_v);
2131  elements[new_e] = new Segment(new_v, vert[1], attr);
2132  vert[1] = new_v;
2133 
2134  CoarseFineTr.embeddings[i] = Embedding(i, 1);
2135  CoarseFineTr.embeddings[new_e] = Embedding(i, 2);
2136  }
2137 
2138  static double seg_children[3*2] = { 0.0,1.0, 0.0,0.5, 0.5,1.0 };
2139  CoarseFineTr.point_matrices.UseExternalData(seg_children, 1, 2, 3);
2140 
2141  GenerateFaces();
2142  } // end of 'if (Dim == 1)'
2143 
2145  sequence++;
2146 
2147  UpdateNodes();
2148 
2149 #ifdef MFEM_DEBUG
2150  CheckElementOrientation(false);
2152 #endif
2153 }
2154 
2156  int nc_limit)
2157 {
2158  if (NURBSext)
2159  {
2160  MFEM_ABORT("ParMesh::NonconformingRefinement: NURBS meshes are not "
2161  "supported. Project the NURBS to Nodes first.");
2162  }
2163 
2164  if (!pncmesh)
2165  {
2166  MFEM_ABORT("Can't convert conforming ParMesh to nonconforming ParMesh "
2167  "(you need to initialize the ParMesh from a nonconforming "
2168  "serial Mesh)");
2169  }
2170 
2171  // NOTE: no check of !refinements.Size(), in parallel we would have to reduce
2172 
2173  // do the refinements
2175  pncmesh->Refine(refinements);
2176 
2177  if (nc_limit > 0)
2178  {
2179  pncmesh->LimitNCLevel(nc_limit);
2180  }
2181 
2182  // create a second mesh containing the finest elements from 'pncmesh'
2183  ParMesh* pmesh2 = new ParMesh(*pncmesh);
2184  pncmesh->OnMeshUpdated(pmesh2);
2185 
2186  attributes.Copy(pmesh2->attributes);
2188 
2189  // now swap the meshes, the second mesh will become the old coarse mesh
2190  // and this mesh will be the new fine mesh
2191  Swap(*pmesh2, false);
2192 
2193  delete pmesh2; // NOTE: old face neighbors destroyed here
2194 
2196 
2198  sequence++;
2199 
2200  if (Nodes) // update/interpolate curved mesh
2201  {
2202  Nodes->FESpace()->Update();
2203  Nodes->Update();
2204  }
2205 }
2206 
2208  double threshold, int nc_limit, int op)
2209 {
2210  const Table &dt = pncmesh->GetDerefinementTable();
2211 
2212  pncmesh->SynchronizeDerefinementData(elem_error, dt);
2213 
2214  Array<int> level_ok;
2215  if (nc_limit > 0)
2216  {
2217  pncmesh->CheckDerefinementNCLevel(dt, level_ok, nc_limit);
2218  }
2219 
2220  Array<int> derefs;
2221  for (int i = 0; i < dt.Size(); i++)
2222  {
2223  if (nc_limit > 0 && !level_ok[i]) { continue; }
2224 
2225  const int* fine = dt.GetRow(i);
2226  int size = dt.RowSize(i);
2227 
2228  double error = 0.0;
2229  for (int j = 0; j < size; j++)
2230  {
2231  MFEM_VERIFY(fine[j] < elem_error.Size(), "");
2232 
2233  double err_fine = elem_error[fine[j]];
2234  switch (op)
2235  {
2236  case 0: error = std::min(error, err_fine); break;
2237  case 1: error += err_fine; break;
2238  case 2: error = std::max(error, err_fine); break;
2239  }
2240  }
2241 
2242  if (error < threshold) { derefs.Append(i); }
2243  }
2244 
2245  long glob_size = ReduceInt(derefs.Size());
2246  if (glob_size)
2247  {
2248  DerefineMesh(derefs);
2249  return true;
2250  }
2251 
2252  return false;
2253 }
2254 
2256 {
2257  if (Conforming())
2258  {
2259  MFEM_ABORT("Load balancing is currently not supported for conforming"
2260  " meshes.");
2261  }
2262 
2264 
2265  pncmesh->Rebalance();
2266 
2267  ParMesh* pmesh2 = new ParMesh(*pncmesh);
2268  pncmesh->OnMeshUpdated(pmesh2);
2269 
2270  attributes.Copy(pmesh2->attributes);
2272 
2273  Swap(*pmesh2, false);
2274  delete pmesh2;
2275 
2277 
2279  sequence++;
2280 
2281  if (Nodes) // redistribute curved mesh
2282  {
2283  Nodes->FESpace()->Update();
2284  Nodes->Update();
2285  }
2286 }
2287 
2288 void ParMesh::RefineGroups(const DSTable &v_to_v, int *middle)
2289 {
2290  int i, attr, newv[3], ind, f_ind, *v;
2291 
2292  int group;
2293  Array<int> group_verts, group_edges, group_faces;
2294 
2295  // To update the groups after a refinement, we observe that:
2296  // - every (new and old) vertex, edge and face belongs to exactly one group
2297  // - the refinement does not create new groups
2298  // - a new vertex appears only as the middle of a refined edge
2299  // - a face can be refined 2, 3 or 4 times producing new edges and faces
2300 
2301  int *I_group_svert, *J_group_svert;
2302  int *I_group_sedge, *J_group_sedge;
2303  int *I_group_sface, *J_group_sface;
2304 
2305  I_group_svert = new int[GetNGroups()+1];
2306  I_group_sedge = new int[GetNGroups()+1];
2307  if (Dim == 3)
2308  {
2309  I_group_sface = new int[GetNGroups()+1];
2310  }
2311  else
2312  {
2313  I_group_sface = NULL;
2314  }
2315 
2316  I_group_svert[0] = I_group_svert[1] = 0;
2317  I_group_sedge[0] = I_group_sedge[1] = 0;
2318  if (Dim == 3)
2319  {
2320  I_group_sface[0] = I_group_sface[1] = 0;
2321  }
2322 
2323  // overestimate the size of the J arrays
2324  if (Dim == 3)
2325  {
2326  J_group_svert = new int[group_svert.Size_of_connections()
2328  J_group_sedge = new int[2*group_sedge.Size_of_connections()
2330  J_group_sface = new int[4*group_sface.Size_of_connections()];
2331  }
2332  else if (Dim == 2)
2333  {
2334  J_group_svert = new int[group_svert.Size_of_connections()
2336  J_group_sedge = new int[2*group_sedge.Size_of_connections()];
2337  J_group_sface = NULL;
2338  }
2339  else
2340  {
2341  J_group_svert = J_group_sedge = J_group_sface = NULL;
2342  }
2343 
2344  for (group = 0; group < GetNGroups()-1; group++)
2345  {
2346  // Get the group shared objects
2347  group_svert.GetRow(group, group_verts);
2348  group_sedge.GetRow(group, group_edges);
2349  group_sface.GetRow(group, group_faces);
2350 
2351  // Check which edges have been refined
2352  for (i = 0; i < group_sedge.RowSize(group); i++)
2353  {
2354  v = shared_edges[group_edges[i]]->GetVertices();
2355  ind = middle[v_to_v(v[0], v[1])];
2356  if (ind != -1)
2357  {
2358  // add a vertex
2359  group_verts.Append(svert_lvert.Append(ind)-1);
2360  // update the edges
2361  attr = shared_edges[group_edges[i]]->GetAttribute();
2362  shared_edges.Append(new Segment(v[1], ind, attr));
2363  group_edges.Append(sedge_ledge.Append(-1)-1);
2364  v[1] = ind;
2365  }
2366  }
2367 
2368  // Check which faces have been refined
2369  for (i = 0; i < group_sface.RowSize(group); i++)
2370  {
2371  v = shared_faces[group_faces[i]]->GetVertices();
2372  ind = middle[v_to_v(v[0], v[1])];
2373  if (ind != -1)
2374  {
2375  attr = shared_faces[group_faces[i]]->GetAttribute();
2376  // add the refinement edge
2377  shared_edges.Append(new Segment(v[2], ind, attr));
2378  group_edges.Append(sedge_ledge.Append(-1)-1);
2379  // add a face
2380  f_ind = group_faces.Size();
2381  shared_faces.Append(new Triangle(v[1], v[2], ind, attr));
2382  group_faces.Append(sface_lface.Append(-1)-1);
2383  newv[0] = v[2]; newv[1] = v[0]; newv[2] = ind;
2384  shared_faces[group_faces[i]]->SetVertices(newv);
2385 
2386  // check if the left face has also been refined
2387  // v = shared_faces[group_faces[i]]->GetVertices();
2388  ind = middle[v_to_v(v[0], v[1])];
2389  if (ind != -1)
2390  {
2391  // add the refinement edge
2392  shared_edges.Append(new Segment(v[2], ind, attr));
2393  group_edges.Append(sedge_ledge.Append(-1)-1);
2394  // add a face
2395  shared_faces.Append(new Triangle(v[1], v[2], ind, attr));
2396  group_faces.Append(sface_lface.Append(-1)-1);
2397  newv[0] = v[2]; newv[1] = v[0]; newv[2] = ind;
2398  shared_faces[group_faces[i]]->SetVertices(newv);
2399  }
2400 
2401  // check if the right face has also been refined
2402  v = shared_faces[group_faces[f_ind]]->GetVertices();
2403  ind = middle[v_to_v(v[0], v[1])];
2404  if (ind != -1)
2405  {
2406  // add the refinement edge
2407  shared_edges.Append(new Segment(v[2], ind, attr));
2408  group_edges.Append(sedge_ledge.Append(-1)-1);
2409  // add a face
2410  shared_faces.Append(new Triangle(v[1], v[2], ind, attr));
2411  group_faces.Append(sface_lface.Append(-1)-1);
2412  newv[0] = v[2]; newv[1] = v[0]; newv[2] = ind;
2413  shared_faces[group_faces[f_ind]]->SetVertices(newv);
2414  }
2415  }
2416  }
2417 
2418  I_group_svert[group+1] = I_group_svert[group] + group_verts.Size();
2419  I_group_sedge[group+1] = I_group_sedge[group] + group_edges.Size();
2420  if (Dim == 3)
2421  {
2422  I_group_sface[group+1] = I_group_sface[group] + group_faces.Size();
2423  }
2424 
2425  int *J;
2426  J = J_group_svert+I_group_svert[group];
2427  for (i = 0; i < group_verts.Size(); i++)
2428  {
2429  J[i] = group_verts[i];
2430  }
2431  J = J_group_sedge+I_group_sedge[group];
2432  for (i = 0; i < group_edges.Size(); i++)
2433  {
2434  J[i] = group_edges[i];
2435  }
2436  if (Dim == 3)
2437  {
2438  J = J_group_sface+I_group_sface[group];
2439  for (i = 0; i < group_faces.Size(); i++)
2440  {
2441  J[i] = group_faces[i];
2442  }
2443  }
2444  }
2445 
2446  // Fix the local numbers of shared edges and faces
2447  {
2448  DSTable new_v_to_v(NumOfVertices);
2449  GetVertexToVertexTable(new_v_to_v);
2450  for (i = 0; i < shared_edges.Size(); i++)
2451  {
2452  v = shared_edges[i]->GetVertices();
2453  sedge_ledge[i] = new_v_to_v(v[0], v[1]);
2454  }
2455  }
2456  if (Dim == 3)
2457  {
2458  STable3D *faces_tbl = GetElementToFaceTable(1);
2459  for (i = 0; i < shared_faces.Size(); i++)
2460  {
2461  v = shared_faces[i]->GetVertices();
2462  sface_lface[i] = (*faces_tbl)(v[0], v[1], v[2]);
2463  }
2464  delete faces_tbl;
2465  }
2466 
2467  group_svert.SetIJ(I_group_svert, J_group_svert);
2468  group_sedge.SetIJ(I_group_sedge, J_group_sedge);
2469  if (Dim == 3)
2470  {
2471  group_sface.SetIJ(I_group_sface, J_group_sface);
2472  }
2473 }
2474 
2476 {
2478 
2479  int oedge = NumOfVertices;
2480 
2481  // call Mesh::QuadUniformRefinement so that it won't update the nodes
2482  {
2483  GridFunction *nodes = Nodes;
2484  Nodes = NULL;
2486  Nodes = nodes;
2487  }
2488 
2489  // update the groups
2490  {
2491  int i, attr, ind, *v;
2492 
2493  int group;
2494  Array<int> sverts, sedges;
2495 
2496  int *I_group_svert, *J_group_svert;
2497  int *I_group_sedge, *J_group_sedge;
2498 
2499  I_group_svert = new int[GetNGroups()+1];
2500  I_group_sedge = new int[GetNGroups()+1];
2501 
2502  I_group_svert[0] = I_group_svert[1] = 0;
2503  I_group_sedge[0] = I_group_sedge[1] = 0;
2504 
2505  // compute the size of the J arrays
2506  J_group_svert = new int[group_svert.Size_of_connections()
2508  J_group_sedge = new int[2*group_sedge.Size_of_connections()];
2509 
2510  for (group = 0; group < GetNGroups()-1; group++)
2511  {
2512  // Get the group shared objects
2513  group_svert.GetRow(group, sverts);
2514  group_sedge.GetRow(group, sedges);
2515 
2516  // Process all the edges
2517  for (i = 0; i < group_sedge.RowSize(group); i++)
2518  {
2519  v = shared_edges[sedges[i]]->GetVertices();
2520  ind = oedge + sedge_ledge[sedges[i]];
2521  // add a vertex
2522  sverts.Append(svert_lvert.Append(ind)-1);
2523  // update the edges
2524  attr = shared_edges[sedges[i]]->GetAttribute();
2525  shared_edges.Append(new Segment(v[1], ind, attr));
2526  sedges.Append(sedge_ledge.Append(-1)-1);
2527  v[1] = ind;
2528  }
2529 
2530  I_group_svert[group+1] = I_group_svert[group] + sverts.Size();
2531  I_group_sedge[group+1] = I_group_sedge[group] + sedges.Size();
2532 
2533  int *J;
2534  J = J_group_svert+I_group_svert[group];
2535  for (i = 0; i < sverts.Size(); i++)
2536  {
2537  J[i] = sverts[i];
2538  }
2539  J = J_group_sedge+I_group_sedge[group];
2540  for (i = 0; i < sedges.Size(); i++)
2541  {
2542  J[i] = sedges[i];
2543  }
2544  }
2545 
2546  // Fix the local numbers of shared edges
2547  DSTable v_to_v(NumOfVertices);
2548  GetVertexToVertexTable(v_to_v);
2549  for (i = 0; i < shared_edges.Size(); i++)
2550  {
2551  v = shared_edges[i]->GetVertices();
2552  sedge_ledge[i] = v_to_v(v[0], v[1]);
2553  }
2554 
2555  group_svert.SetIJ(I_group_svert, J_group_svert);
2556  group_sedge.SetIJ(I_group_sedge, J_group_sedge);
2557  }
2558 
2559  UpdateNodes();
2560 }
2561 
2563 {
2565 
2566  int oedge = NumOfVertices;
2567  int oface = oedge + NumOfEdges;
2568 
2569  DSTable v_to_v(NumOfVertices);
2570  GetVertexToVertexTable(v_to_v);
2571  STable3D *faces_tbl = GetFacesTable();
2572 
2573  // call Mesh::HexUniformRefinement so that it won't update the nodes
2574  {
2575  GridFunction *nodes = Nodes;
2576  Nodes = NULL;
2578  Nodes = nodes;
2579  }
2580 
2581  // update the groups
2582  {
2583  int i, attr, newv[4], ind, m[5];
2584  Array<int> v;
2585 
2586  int group;
2587  Array<int> group_verts, group_edges, group_faces;
2588 
2589  int *I_group_svert, *J_group_svert;
2590  int *I_group_sedge, *J_group_sedge;
2591  int *I_group_sface, *J_group_sface;
2592 
2593  I_group_svert = new int[GetNGroups()+1];
2594  I_group_sedge = new int[GetNGroups()+1];
2595  I_group_sface = new int[GetNGroups()+1];
2596 
2597  I_group_svert[0] = I_group_svert[1] = 0;
2598  I_group_sedge[0] = I_group_sedge[1] = 0;
2599  I_group_sface[0] = I_group_sface[1] = 0;
2600 
2601  // compute the size of the J arrays
2602  J_group_svert = new int[group_svert.Size_of_connections()
2605  J_group_sedge = new int[2*group_sedge.Size_of_connections()
2607  J_group_sface = new int[4*group_sface.Size_of_connections()];
2608 
2609  for (group = 0; group < GetNGroups()-1; group++)
2610  {
2611  // Get the group shared objects
2612  group_svert.GetRow(group, group_verts);
2613  group_sedge.GetRow(group, group_edges);
2614  group_sface.GetRow(group, group_faces);
2615 
2616  // Process the edges that have been refined
2617  for (i = 0; i < group_sedge.RowSize(group); i++)
2618  {
2619  shared_edges[group_edges[i]]->GetVertices(v);
2620  ind = oedge + v_to_v(v[0], v[1]);
2621  // add a vertex
2622  group_verts.Append(svert_lvert.Append(ind)-1);
2623  // update the edges
2624  attr = shared_edges[group_edges[i]]->GetAttribute();
2625  shared_edges.Append(new Segment(v[1], ind, attr));
2626  group_edges.Append(sedge_ledge.Append(-1)-1);
2627  newv[0] = v[0]; newv[1] = ind;
2628  shared_edges[group_edges[i]]->SetVertices(newv);
2629  }
2630 
2631  // Process the faces that have been refined
2632  for (i = 0; i < group_sface.RowSize(group); i++)
2633  {
2634  shared_faces[group_faces[i]]->GetVertices(v);
2635  m[0] = oface+(*faces_tbl)(v[0], v[1], v[2], v[3]);
2636  // add a vertex
2637  group_verts.Append(svert_lvert.Append(m[0])-1);
2638  // add the refinement edges
2639  attr = shared_faces[group_faces[i]]->GetAttribute();
2640  m[1] = oedge + v_to_v(v[0], v[1]);
2641  m[2] = oedge + v_to_v(v[1], v[2]);
2642  m[3] = oedge + v_to_v(v[2], v[3]);
2643  m[4] = oedge + v_to_v(v[3], v[0]);
2644  shared_edges.Append(new Segment(m[1], m[0], attr));
2645  group_edges.Append(sedge_ledge.Append(-1)-1);
2646  shared_edges.Append(new Segment(m[2], m[0], attr));
2647  group_edges.Append(sedge_ledge.Append(-1)-1);
2648  shared_edges.Append(new Segment(m[3], m[0], attr));
2649  group_edges.Append(sedge_ledge.Append(-1)-1);
2650  shared_edges.Append(new Segment(m[4], m[0], attr));
2651  group_edges.Append(sedge_ledge.Append(-1)-1);
2652  // update faces
2653  newv[0] = v[0]; newv[1] = m[1]; newv[2] = m[0]; newv[3] = m[4];
2654  shared_faces[group_faces[i]]->SetVertices(newv);
2655  shared_faces.Append(new Quadrilateral(m[1],v[1],m[2],m[0],attr));
2656  group_faces.Append(sface_lface.Append(-1)-1);
2657  shared_faces.Append(new Quadrilateral(m[0],m[2],v[2],m[3],attr));
2658  group_faces.Append(sface_lface.Append(-1)-1);
2659  shared_faces.Append(new Quadrilateral(m[4],m[0],m[3],v[3],attr));
2660  group_faces.Append(sface_lface.Append(-1)-1);
2661  }
2662 
2663  I_group_svert[group+1] = I_group_svert[group] + group_verts.Size();
2664  I_group_sedge[group+1] = I_group_sedge[group] + group_edges.Size();
2665  I_group_sface[group+1] = I_group_sface[group] + group_faces.Size();
2666 
2667  int *J;
2668  J = J_group_svert+I_group_svert[group];
2669  for (i = 0; i < group_verts.Size(); i++)
2670  {
2671  J[i] = group_verts[i];
2672  }
2673  J = J_group_sedge+I_group_sedge[group];
2674  for (i = 0; i < group_edges.Size(); i++)
2675  {
2676  J[i] = group_edges[i];
2677  }
2678  J = J_group_sface+I_group_sface[group];
2679  for (i = 0; i < group_faces.Size(); i++)
2680  {
2681  J[i] = group_faces[i];
2682  }
2683  }
2684 
2685  // Fix the local numbers of shared edges and faces
2686  DSTable new_v_to_v(NumOfVertices);
2687  GetVertexToVertexTable(new_v_to_v);
2688  for (i = 0; i < shared_edges.Size(); i++)
2689  {
2690  shared_edges[i]->GetVertices(v);
2691  sedge_ledge[i] = new_v_to_v(v[0], v[1]);
2692  }
2693 
2694  delete faces_tbl;
2695  faces_tbl = GetFacesTable();
2696  for (i = 0; i < shared_faces.Size(); i++)
2697  {
2698  shared_faces[i]->GetVertices(v);
2699  sface_lface[i] = (*faces_tbl)(v[0], v[1], v[2], v[3]);
2700  }
2701  delete faces_tbl;
2702 
2703  group_svert.SetIJ(I_group_svert, J_group_svert);
2704  group_sedge.SetIJ(I_group_sedge, J_group_sedge);
2705  group_sface.SetIJ(I_group_sface, J_group_sface);
2706  }
2707 
2708  UpdateNodes();
2709 }
2710 
2712 {
2713  if (MyRank == 0)
2714  {
2715  cout << "\nParMesh::NURBSUniformRefinement : Not supported yet!\n";
2716  }
2717 }
2718 
2719 void ParMesh::PrintXG(std::ostream &out) const
2720 {
2721  MFEM_ASSERT(Dim == spaceDim, "2D manifolds not supported");
2722  if (Dim == 3 && meshgen == 1)
2723  {
2724  int i, j, nv;
2725  const int *ind;
2726 
2727  out << "NETGEN_Neutral_Format\n";
2728  // print the vertices
2729  out << NumOfVertices << '\n';
2730  for (i = 0; i < NumOfVertices; i++)
2731  {
2732  for (j = 0; j < Dim; j++)
2733  {
2734  out << " " << vertices[i](j);
2735  }
2736  out << '\n';
2737  }
2738 
2739  // print the elements
2740  out << NumOfElements << '\n';
2741  for (i = 0; i < NumOfElements; i++)
2742  {
2743  nv = elements[i]->GetNVertices();
2744  ind = elements[i]->GetVertices();
2745  out << elements[i]->GetAttribute();
2746  for (j = 0; j < nv; j++)
2747  {
2748  out << " " << ind[j]+1;
2749  }
2750  out << '\n';
2751  }
2752 
2753  // print the boundary + shared faces information
2754  out << NumOfBdrElements + shared_faces.Size() << '\n';
2755  // boundary
2756  for (i = 0; i < NumOfBdrElements; i++)
2757  {
2758  nv = boundary[i]->GetNVertices();
2759  ind = boundary[i]->GetVertices();
2760  out << boundary[i]->GetAttribute();
2761  for (j = 0; j < nv; j++)
2762  {
2763  out << " " << ind[j]+1;
2764  }
2765  out << '\n';
2766  }
2767  // shared faces
2768  for (i = 0; i < shared_faces.Size(); i++)
2769  {
2770  nv = shared_faces[i]->GetNVertices();
2771  ind = shared_faces[i]->GetVertices();
2772  out << shared_faces[i]->GetAttribute();
2773  for (j = 0; j < nv; j++)
2774  {
2775  out << " " << ind[j]+1;
2776  }
2777  out << '\n';
2778  }
2779  }
2780 
2781  if (Dim == 3 && meshgen == 2)
2782  {
2783  int i, j, nv;
2784  const int *ind;
2785 
2786  out << "TrueGrid\n"
2787  << "1 " << NumOfVertices << " " << NumOfElements << " 0 0 0 0 0 0 0\n"
2788  << "0 0 0 1 0 0 0 0 0 0 0\n"
2789  << "0 0 " << NumOfBdrElements+shared_faces.Size()
2790  << " 0 0 0 0 0 0 0 0 0 0 0 0 0\n"
2791  << "0.0 0.0 0.0 0 0 0.0 0.0 0 0.0\n"
2792  << "0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0\n";
2793 
2794  // print the vertices
2795  for (i = 0; i < NumOfVertices; i++)
2796  out << i+1 << " 0.0 " << vertices[i](0) << " " << vertices[i](1)
2797  << " " << vertices[i](2) << " 0.0\n";
2798 
2799  // print the elements
2800  for (i = 0; i < NumOfElements; i++)
2801  {
2802  nv = elements[i]->GetNVertices();
2803  ind = elements[i]->GetVertices();
2804  out << i+1 << " " << elements[i]->GetAttribute();
2805  for (j = 0; j < nv; j++)
2806  {
2807  out << " " << ind[j]+1;
2808  }
2809  out << '\n';
2810  }
2811 
2812  // print the boundary information
2813  for (i = 0; i < NumOfBdrElements; i++)
2814  {
2815  nv = boundary[i]->GetNVertices();
2816  ind = boundary[i]->GetVertices();
2817  out << boundary[i]->GetAttribute();
2818  for (j = 0; j < nv; j++)
2819  {
2820  out << " " << ind[j]+1;
2821  }
2822  out << " 1.0 1.0 1.0 1.0\n";
2823  }
2824 
2825  // print the shared faces information
2826  for (i = 0; i < shared_faces.Size(); i++)
2827  {
2828  nv = shared_faces[i]->GetNVertices();
2829  ind = shared_faces[i]->GetVertices();
2830  out << shared_faces[i]->GetAttribute();
2831  for (j = 0; j < nv; j++)
2832  {
2833  out << " " << ind[j]+1;
2834  }
2835  out << " 1.0 1.0 1.0 1.0\n";
2836  }
2837  }
2838 
2839  if (Dim == 2)
2840  {
2841  int i, j, attr;
2842  Array<int> v;
2843 
2844  out << "areamesh2\n\n";
2845 
2846  // print the boundary + shared edges information
2847  out << NumOfBdrElements + shared_edges.Size() << '\n';
2848  // boundary
2849  for (i = 0; i < NumOfBdrElements; i++)
2850  {
2851  attr = boundary[i]->GetAttribute();
2852  boundary[i]->GetVertices(v);
2853  out << attr << " ";
2854  for (j = 0; j < v.Size(); j++)
2855  {
2856  out << v[j] + 1 << " ";
2857  }
2858  out << '\n';
2859  }
2860  // shared edges
2861  for (i = 0; i < shared_edges.Size(); i++)
2862  {
2863  attr = shared_edges[i]->GetAttribute();
2864  shared_edges[i]->GetVertices(v);
2865  out << attr << " ";
2866  for (j = 0; j < v.Size(); j++)
2867  {
2868  out << v[j] + 1 << " ";
2869  }
2870  out << '\n';
2871  }
2872 
2873  // print the elements
2874  out << NumOfElements << '\n';
2875  for (i = 0; i < NumOfElements; i++)
2876  {
2877  attr = elements[i]->GetAttribute();
2878  elements[i]->GetVertices(v);
2879 
2880  out << attr << " ";
2881  if ((j = GetElementType(i)) == Element::TRIANGLE)
2882  {
2883  out << 3 << " ";
2884  }
2885  else if (j == Element::QUADRILATERAL)
2886  {
2887  out << 4 << " ";
2888  }
2889  else if (j == Element::SEGMENT)
2890  {
2891  out << 2 << " ";
2892  }
2893  for (j = 0; j < v.Size(); j++)
2894  {
2895  out << v[j] + 1 << " ";
2896  }
2897  out << '\n';
2898  }
2899 
2900  // print the vertices
2901  out << NumOfVertices << '\n';
2902  for (i = 0; i < NumOfVertices; i++)
2903  {
2904  for (j = 0; j < Dim; j++)
2905  {
2906  out << vertices[i](j) << " ";
2907  }
2908  out << '\n';
2909  }
2910  }
2911  out.flush();
2912 }
2913 
2915 {
2916  // In 2D, this is a workaround for a CPU boundary rendering artifact. We need
2917  // to skip a shared master edge if one of its slaves has the same rank.
2918 
2919  const NCMesh::NCList &list = pncmesh->GetEdgeList();
2920  for (int i = master.slaves_begin; i < master.slaves_end; i++)
2921  {
2922  if (!pncmesh->IsGhost(1, list.slaves[i].index)) { return true; }
2923  }
2924  return false;
2925 }
2926 
2927 void ParMesh::Print(std::ostream &out) const
2928 {
2929  bool print_shared = true;
2930  int i, j, shared_bdr_attr;
2931  Array<int> nc_shared_faces;
2932 
2933  const Array<int>* s2l_face;
2934  if (!pncmesh)
2935  {
2936  s2l_face = ((Dim == 1) ? &svert_lvert :
2937  ((Dim == 2) ? &sedge_ledge : &sface_lface));
2938  }
2939  else
2940  {
2941  s2l_face = &nc_shared_faces;
2942  if (Dim >= 2)
2943  {
2944  // get a list of all shared non-ghost faces
2945  const NCMesh::NCList& sfaces =
2946  (Dim == 3) ? pncmesh->GetSharedFaces() : pncmesh->GetSharedEdges();
2947  const int nfaces = GetNumFaces();
2948  for (unsigned i = 0; i < sfaces.conforming.size(); i++)
2949  {
2950  int index = sfaces.conforming[i].index;
2951  if (index < nfaces) { nc_shared_faces.Append(index); }
2952  }
2953  for (unsigned i = 0; i < sfaces.masters.size(); i++)
2954  {
2955  if (Dim == 2 && WantSkipSharedMaster(sfaces.masters[i])) { continue; }
2956  int index = sfaces.masters[i].index;
2957  if (index < nfaces) { nc_shared_faces.Append(index); }
2958  }
2959  for (unsigned i = 0; i < sfaces.slaves.size(); i++)
2960  {
2961  int index = sfaces.slaves[i].index;
2962  if (index < nfaces) { nc_shared_faces.Append(index); }
2963  }
2964  }
2965  }
2966 
2967  if (NURBSext)
2968  {
2969  Mesh::Print(out); // does not print shared boundary
2970  return;
2971  }
2972 
2973  out << "MFEM mesh v1.0\n";
2974 
2975  // optional
2976  out <<
2977  "\n#\n# MFEM Geometry Types (see mesh/geom.hpp):\n#\n"
2978  "# POINT = 0\n"
2979  "# SEGMENT = 1\n"
2980  "# TRIANGLE = 2\n"
2981  "# SQUARE = 3\n"
2982  "# TETRAHEDRON = 4\n"
2983  "# CUBE = 5\n"
2984  "#\n";
2985 
2986  out << "\ndimension\n" << Dim
2987  << "\n\nelements\n" << NumOfElements << '\n';
2988  for (i = 0; i < NumOfElements; i++)
2989  {
2990  PrintElement(elements[i], out);
2991  }
2992 
2993  int num_bdr_elems = NumOfBdrElements;
2994  if (print_shared && Dim > 1)
2995  {
2996  num_bdr_elems += s2l_face->Size();
2997  }
2998  out << "\nboundary\n" << num_bdr_elems << '\n';
2999  for (i = 0; i < NumOfBdrElements; i++)
3000  {
3001  PrintElement(boundary[i], out);
3002  }
3003 
3004  if (print_shared && Dim > 1)
3005  {
3006  if (bdr_attributes.Size())
3007  {
3008  shared_bdr_attr = bdr_attributes.Max() + MyRank + 1;
3009  }
3010  else
3011  {
3012  shared_bdr_attr = MyRank + 1;
3013  }
3014  for (i = 0; i < s2l_face->Size(); i++)
3015  {
3016  // Modify the attributes of the faces (not used otherwise?)
3017  faces[(*s2l_face)[i]]->SetAttribute(shared_bdr_attr);
3018  PrintElement(faces[(*s2l_face)[i]], out);
3019  }
3020  }
3021  out << "\nvertices\n" << NumOfVertices << '\n';
3022  if (Nodes == NULL)
3023  {
3024  out << spaceDim << '\n';
3025  for (i = 0; i < NumOfVertices; i++)
3026  {
3027  out << vertices[i](0);
3028  for (j = 1; j < spaceDim; j++)
3029  {
3030  out << ' ' << vertices[i](j);
3031  }
3032  out << '\n';
3033  }
3034  out.flush();
3035  }
3036  else
3037  {
3038  out << "\nnodes\n";
3039  Nodes->Save(out);
3040  }
3041 }
3042 
3043 static void dump_element(const Element* elem, Array<int> &data)
3044 {
3045  data.Append(elem->GetGeometryType());
3046 
3047  int nv = elem->GetNVertices();
3048  const int *v = elem->GetVertices();
3049  for (int i = 0; i < nv; i++)
3050  {
3051  data.Append(v[i]);
3052  }
3053 }
3054 
3055 void ParMesh::PrintAsOne(std::ostream &out)
3056 {
3057  int i, j, k, p, nv_ne[2], &nv = nv_ne[0], &ne = nv_ne[1], vc;
3058  const int *v;
3059  MPI_Status status;
3060  Array<double> vert;
3061  Array<int> ints;
3062 
3063  if (MyRank == 0)
3064  {
3065  out << "MFEM mesh v1.0\n";
3066 
3067  // optional
3068  out <<
3069  "\n#\n# MFEM Geometry Types (see mesh/geom.hpp):\n#\n"
3070  "# POINT = 0\n"
3071  "# SEGMENT = 1\n"
3072  "# TRIANGLE = 2\n"
3073  "# SQUARE = 3\n"
3074  "# TETRAHEDRON = 4\n"
3075  "# CUBE = 5\n"
3076  "#\n";
3077 
3078  out << "\ndimension\n" << Dim;
3079  }
3080 
3081  nv = NumOfElements;
3082  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3083  if (MyRank == 0)
3084  {
3085  out << "\n\nelements\n" << ne << '\n';
3086  for (i = 0; i < NumOfElements; i++)
3087  {
3088  // processor number + 1 as attribute and geometry type
3089  out << 1 << ' ' << elements[i]->GetGeometryType();
3090  // vertices
3091  nv = elements[i]->GetNVertices();
3092  v = elements[i]->GetVertices();
3093  for (j = 0; j < nv; j++)
3094  {
3095  out << ' ' << v[j];
3096  }
3097  out << '\n';
3098  }
3099  vc = NumOfVertices;
3100  for (p = 1; p < NRanks; p++)
3101  {
3102  MPI_Recv(nv_ne, 2, MPI_INT, p, 444, MyComm, &status);
3103  ints.SetSize(ne);
3104  if (ne)
3105  {
3106  MPI_Recv(&ints[0], ne, MPI_INT, p, 445, MyComm, &status);
3107  }
3108  for (i = 0; i < ne; )
3109  {
3110  // processor number + 1 as attribute and geometry type
3111  out << p+1 << ' ' << ints[i];
3112  // vertices
3113  k = Geometries.GetVertices(ints[i++])->GetNPoints();
3114  for (j = 0; j < k; j++)
3115  {
3116  out << ' ' << vc + ints[i++];
3117  }
3118  out << '\n';
3119  }
3120  vc += nv;
3121  }
3122  }
3123  else
3124  {
3125  // for each element send its geometry type and its vertices
3126  ne = 0;
3127  for (i = 0; i < NumOfElements; i++)
3128  {
3129  ne += 1 + elements[i]->GetNVertices();
3130  }
3131  nv = NumOfVertices;
3132  MPI_Send(nv_ne, 2, MPI_INT, 0, 444, MyComm);
3133 
3134  ints.Reserve(ne);
3135  ints.SetSize(0);
3136  for (i = 0; i < NumOfElements; i++)
3137  {
3138  dump_element(elements[i], ints);
3139  }
3140  MFEM_ASSERT(ints.Size() == ne, "");
3141  if (ne)
3142  {
3143  MPI_Send(&ints[0], ne, MPI_INT, 0, 445, MyComm);
3144  }
3145  }
3146 
3147  // boundary + shared boundary
3148  ne = NumOfBdrElements;
3149  if (!pncmesh)
3150  {
3151  ne += ((Dim == 2) ? shared_edges : shared_faces).Size();
3152  }
3153  else if (Dim > 1)
3154  {
3155  const NCMesh::NCList &list = pncmesh->GetSharedList(Dim - 1);
3156  ne += list.conforming.size() + list.masters.size() + list.slaves.size();
3157  }
3158  ints.Reserve(ne * (1 + 2*(Dim-1))); // just an upper bound
3159  ints.SetSize(0);
3160 
3161  // for each boundary and shared boundary element send its geometry type
3162  // and its vertices
3163  ne = 0;
3164  for (i = j = 0; i < NumOfBdrElements; i++)
3165  {
3166  dump_element(boundary[i], ints); ne++;
3167  }
3168  if (!pncmesh)
3169  {
3170  Array<Element*> &shared = (Dim == 2) ? shared_edges : shared_faces;
3171  for (i = 0; i < shared.Size(); i++)
3172  {
3173  dump_element(shared[i], ints); ne++;
3174  }
3175  }
3176  else if (Dim > 1)
3177  {
3178  const NCMesh::NCList &list = pncmesh->GetSharedList(Dim - 1);
3179  const int nfaces = GetNumFaces();
3180  for (i = 0; i < (int) list.conforming.size(); i++)
3181  {
3182  int index = list.conforming[i].index;
3183  if (index < nfaces) { dump_element(faces[index], ints); ne++; }
3184  }
3185  for (i = 0; i < (int) list.masters.size(); i++)
3186  {
3187  int index = list.masters[i].index;
3188  if (index < nfaces) { dump_element(faces[index], ints); ne++; }
3189  }
3190  for (i = 0; i < (int) list.slaves.size(); i++)
3191  {
3192  int index = list.slaves[i].index;
3193  if (index < nfaces) { dump_element(faces[index], ints); ne++; }
3194  }
3195  }
3196 
3197  MPI_Reduce(&ne, &k, 1, MPI_INT, MPI_SUM, 0, MyComm);
3198  if (MyRank == 0)
3199  {
3200  out << "\nboundary\n" << k << '\n';
3201  vc = 0;
3202  for (p = 0; p < NRanks; p++)
3203  {
3204  if (p)
3205  {
3206  MPI_Recv(nv_ne, 2, MPI_INT, p, 446, MyComm, &status);
3207  ints.SetSize(ne);
3208  if (ne)
3209  {
3210  MPI_Recv(ints.GetData(), ne, MPI_INT, p, 447, MyComm, &status);
3211  }
3212  }
3213  else
3214  {
3215  ne = ints.Size();
3216  nv = NumOfVertices;
3217  }
3218  for (i = 0; i < ne; )
3219  {
3220  // processor number + 1 as bdr. attr. and bdr. geometry type
3221  out << p+1 << ' ' << ints[i];
3222  k = Geometries.GetVertices(ints[i++])->GetNPoints();
3223  // vertices
3224  for (j = 0; j < k; j++)
3225  {
3226  out << ' ' << vc + ints[i++];
3227  }
3228  out << '\n';
3229  }
3230  vc += nv;
3231  }
3232  }
3233  else
3234  {
3235  nv = NumOfVertices;
3236  ne = ints.Size();
3237  MPI_Send(nv_ne, 2, MPI_INT, 0, 446, MyComm);
3238  if (ne)
3239  {
3240  MPI_Send(ints.GetData(), ne, MPI_INT, 0, 447, MyComm);
3241  }
3242  }
3243 
3244  // vertices / nodes
3245  MPI_Reduce(&NumOfVertices, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3246  if (MyRank == 0)
3247  {
3248  out << "\nvertices\n" << nv << '\n';
3249  }
3250  if (Nodes == NULL)
3251  {
3252  if (MyRank == 0)
3253  {
3254  out << spaceDim << '\n';
3255  for (i = 0; i < NumOfVertices; i++)
3256  {
3257  out << vertices[i](0);
3258  for (j = 1; j < spaceDim; j++)
3259  {
3260  out << ' ' << vertices[i](j);
3261  }
3262  out << '\n';
3263  }
3264  for (p = 1; p < NRanks; p++)
3265  {
3266  MPI_Recv(&nv, 1, MPI_INT, p, 448, MyComm, &status);
3267  vert.SetSize(nv*spaceDim);
3268  if (nv)
3269  {
3270  MPI_Recv(&vert[0], nv*spaceDim, MPI_DOUBLE, p, 449, MyComm, &status);
3271  }
3272  for (i = 0; i < nv; i++)
3273  {
3274  out << vert[i*spaceDim];
3275  for (j = 1; j < spaceDim; j++)
3276  {
3277  out << ' ' << vert[i*spaceDim+j];
3278  }
3279  out << '\n';
3280  }
3281  }
3282  out.flush();
3283  }
3284  else
3285  {
3286  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 448, MyComm);
3288  for (i = 0; i < NumOfVertices; i++)
3289  {
3290  for (j = 0; j < spaceDim; j++)
3291  {
3292  vert[i*spaceDim+j] = vertices[i](j);
3293  }
3294  }
3295  if (NumOfVertices)
3296  {
3297  MPI_Send(&vert[0], NumOfVertices*spaceDim, MPI_DOUBLE, 0, 449, MyComm);
3298  }
3299  }
3300  }
3301  else
3302  {
3303  if (MyRank == 0)
3304  {
3305  out << "\nnodes\n";
3306  }
3307  ParGridFunction *pnodes = dynamic_cast<ParGridFunction *>(Nodes);
3308  if (pnodes)
3309  {
3310  pnodes->SaveAsOne(out);
3311  }
3312  else
3313  {
3314  ParFiniteElementSpace *pfes =
3315  dynamic_cast<ParFiniteElementSpace *>(Nodes->FESpace());
3316  if (pfes)
3317  {
3318  // create a wrapper ParGridFunction
3319  ParGridFunction ParNodes(pfes, Nodes);
3320  ParNodes.SaveAsOne(out);
3321  }
3322  else
3323  {
3324  mfem_error("ParMesh::PrintAsOne : Nodes have no parallel info!");
3325  }
3326  }
3327  }
3328 }
3329 
3330 void ParMesh::PrintAsOneXG(std::ostream &out)
3331 {
3332  MFEM_ASSERT(Dim == spaceDim, "2D Manifolds not supported.");
3333  if (Dim == 3 && meshgen == 1)
3334  {
3335  int i, j, k, nv, ne, p;
3336  const int *ind, *v;
3337  MPI_Status status;
3338  Array<double> vert;
3339  Array<int> ints;
3340 
3341  if (MyRank == 0)
3342  {
3343  out << "NETGEN_Neutral_Format\n";
3344  // print the vertices
3345  ne = NumOfVertices;
3346  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3347  out << nv << '\n';
3348  for (i = 0; i < NumOfVertices; i++)
3349  {
3350  for (j = 0; j < Dim; j++)
3351  {
3352  out << " " << vertices[i](j);
3353  }
3354  out << '\n';
3355  }
3356  for (p = 1; p < NRanks; p++)
3357  {
3358  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3359  vert.SetSize(Dim*nv);
3360  MPI_Recv(&vert[0], Dim*nv, MPI_DOUBLE, p, 445, MyComm, &status);
3361  for (i = 0; i < nv; i++)
3362  {
3363  for (j = 0; j < Dim; j++)
3364  {
3365  out << " " << vert[Dim*i+j];
3366  }
3367  out << '\n';
3368  }
3369  }
3370 
3371  // print the elements
3372  nv = NumOfElements;
3373  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3374  out << ne << '\n';
3375  for (i = 0; i < NumOfElements; i++)
3376  {
3377  nv = elements[i]->GetNVertices();
3378  ind = elements[i]->GetVertices();
3379  out << 1;
3380  for (j = 0; j < nv; j++)
3381  {
3382  out << " " << ind[j]+1;
3383  }
3384  out << '\n';
3385  }
3386  k = NumOfVertices;
3387  for (p = 1; p < NRanks; p++)
3388  {
3389  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3390  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3391  ints.SetSize(4*ne);
3392  MPI_Recv(&ints[0], 4*ne, MPI_INT, p, 447, MyComm, &status);
3393  for (i = 0; i < ne; i++)
3394  {
3395  out << p+1;
3396  for (j = 0; j < 4; j++)
3397  {
3398  out << " " << k+ints[i*4+j]+1;
3399  }
3400  out << '\n';
3401  }
3402  k += nv;
3403  }
3404  // print the boundary + shared faces information
3405  nv = NumOfBdrElements + shared_faces.Size();
3406  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3407  out << ne << '\n';
3408  // boundary
3409  for (i = 0; i < NumOfBdrElements; i++)
3410  {
3411  nv = boundary[i]->GetNVertices();
3412  ind = boundary[i]->GetVertices();
3413  out << 1;
3414  for (j = 0; j < nv; j++)
3415  {
3416  out << " " << ind[j]+1;
3417  }
3418  out << '\n';
3419  }
3420  // shared faces
3421  for (i = 0; i < shared_faces.Size(); i++)
3422  {
3423  nv = shared_faces[i]->GetNVertices();
3424  ind = shared_faces[i]->GetVertices();
3425  out << 1;
3426  for (j = 0; j < nv; j++)
3427  {
3428  out << " " << ind[j]+1;
3429  }
3430  out << '\n';
3431  }
3432  k = NumOfVertices;
3433  for (p = 1; p < NRanks; p++)
3434  {
3435  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3436  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3437  ints.SetSize(3*ne);
3438  MPI_Recv(&ints[0], 3*ne, MPI_INT, p, 447, MyComm, &status);
3439  for (i = 0; i < ne; i++)
3440  {
3441  out << p+1;
3442  for (j = 0; j < 3; j++)
3443  {
3444  out << " " << k+ints[i*3+j]+1;
3445  }
3446  out << '\n';
3447  }
3448  k += nv;
3449  }
3450  }
3451  else
3452  {
3453  ne = NumOfVertices;
3454  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3455  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3456  vert.SetSize(Dim*NumOfVertices);
3457  for (i = 0; i < NumOfVertices; i++)
3458  for (j = 0; j < Dim; j++)
3459  {
3460  vert[Dim*i+j] = vertices[i](j);
3461  }
3462  MPI_Send(&vert[0], Dim*NumOfVertices, MPI_DOUBLE,
3463  0, 445, MyComm);
3464  // elements
3465  ne = NumOfElements;
3466  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3467  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3468  MPI_Send(&NumOfElements, 1, MPI_INT, 0, 446, MyComm);
3469  ints.SetSize(NumOfElements*4);
3470  for (i = 0; i < NumOfElements; i++)
3471  {
3472  v = elements[i]->GetVertices();
3473  for (j = 0; j < 4; j++)
3474  {
3475  ints[4*i+j] = v[j];
3476  }
3477  }
3478  MPI_Send(&ints[0], 4*NumOfElements, MPI_INT, 0, 447, MyComm);
3479  // boundary + shared faces
3480  nv = NumOfBdrElements + shared_faces.Size();
3481  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3482  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3483  ne = NumOfBdrElements + shared_faces.Size();
3484  MPI_Send(&ne, 1, MPI_INT, 0, 446, MyComm);
3485  ints.SetSize(3*ne);
3486  for (i = 0; i < NumOfBdrElements; i++)
3487  {
3488  v = boundary[i]->GetVertices();
3489  for (j = 0; j < 3; j++)
3490  {
3491  ints[3*i+j] = v[j];
3492  }
3493  }
3494  for ( ; i < ne; i++)
3495  {
3496  v = shared_faces[i-NumOfBdrElements]->GetVertices();
3497  for (j = 0; j < 3; j++)
3498  {
3499  ints[3*i+j] = v[j];
3500  }
3501  }
3502  MPI_Send(&ints[0], 3*ne, MPI_INT, 0, 447, MyComm);
3503  }
3504  }
3505 
3506  if (Dim == 3 && meshgen == 2)
3507  {
3508  int i, j, k, nv, ne, p;
3509  const int *ind, *v;
3510  MPI_Status status;
3511  Array<double> vert;
3512  Array<int> ints;
3513 
3514  int TG_nv, TG_ne, TG_nbe;
3515 
3516  if (MyRank == 0)
3517  {
3518  MPI_Reduce(&NumOfVertices, &TG_nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3519  MPI_Reduce(&NumOfElements, &TG_ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3520  nv = NumOfBdrElements + shared_faces.Size();
3521  MPI_Reduce(&nv, &TG_nbe, 1, MPI_INT, MPI_SUM, 0, MyComm);
3522 
3523  out << "TrueGrid\n"
3524  << "1 " << TG_nv << " " << TG_ne << " 0 0 0 0 0 0 0\n"
3525  << "0 0 0 1 0 0 0 0 0 0 0\n"
3526  << "0 0 " << TG_nbe << " 0 0 0 0 0 0 0 0 0 0 0 0 0\n"
3527  << "0.0 0.0 0.0 0 0 0.0 0.0 0 0.0\n"
3528  << "0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0\n";
3529 
3530  // print the vertices
3531  nv = TG_nv;
3532  for (i = 0; i < NumOfVertices; i++)
3533  out << i+1 << " 0.0 " << vertices[i](0) << " " << vertices[i](1)
3534  << " " << vertices[i](2) << " 0.0\n";
3535  for (p = 1; p < NRanks; p++)
3536  {
3537  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3538  vert.SetSize(Dim*nv);
3539  MPI_Recv(&vert[0], Dim*nv, MPI_DOUBLE, p, 445, MyComm, &status);
3540  for (i = 0; i < nv; i++)
3541  out << i+1 << " 0.0 " << vert[Dim*i] << " " << vert[Dim*i+1]
3542  << " " << vert[Dim*i+2] << " 0.0\n";
3543  }
3544 
3545  // print the elements
3546  ne = TG_ne;
3547  for (i = 0; i < NumOfElements; i++)
3548  {
3549  nv = elements[i]->GetNVertices();
3550  ind = elements[i]->GetVertices();
3551  out << i+1 << " " << 1;
3552  for (j = 0; j < nv; j++)
3553  {
3554  out << " " << ind[j]+1;
3555  }
3556  out << '\n';
3557  }
3558  k = NumOfVertices;
3559  for (p = 1; p < NRanks; p++)
3560  {
3561  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3562  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3563  ints.SetSize(8*ne);
3564  MPI_Recv(&ints[0], 8*ne, MPI_INT, p, 447, MyComm, &status);
3565  for (i = 0; i < ne; i++)
3566  {
3567  out << i+1 << " " << p+1;
3568  for (j = 0; j < 8; j++)
3569  {
3570  out << " " << k+ints[i*8+j]+1;
3571  }
3572  out << '\n';
3573  }
3574  k += nv;
3575  }
3576 
3577  // print the boundary + shared faces information
3578  ne = TG_nbe;
3579  // boundary
3580  for (i = 0; i < NumOfBdrElements; i++)
3581  {
3582  nv = boundary[i]->GetNVertices();
3583  ind = boundary[i]->GetVertices();
3584  out << 1;
3585  for (j = 0; j < nv; j++)
3586  {
3587  out << " " << ind[j]+1;
3588  }
3589  out << " 1.0 1.0 1.0 1.0\n";
3590  }
3591  // shared faces
3592  for (i = 0; i < shared_faces.Size(); i++)
3593  {
3594  nv = shared_faces[i]->GetNVertices();
3595  ind = shared_faces[i]->GetVertices();
3596  out << 1;
3597  for (j = 0; j < nv; j++)
3598  {
3599  out << " " << ind[j]+1;
3600  }
3601  out << " 1.0 1.0 1.0 1.0\n";
3602  }
3603  k = NumOfVertices;
3604  for (p = 1; p < NRanks; p++)
3605  {
3606  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3607  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3608  ints.SetSize(4*ne);
3609  MPI_Recv(&ints[0], 4*ne, MPI_INT, p, 447, MyComm, &status);
3610  for (i = 0; i < ne; i++)
3611  {
3612  out << p+1;
3613  for (j = 0; j < 4; j++)
3614  {
3615  out << " " << k+ints[i*4+j]+1;
3616  }
3617  out << " 1.0 1.0 1.0 1.0\n";
3618  }
3619  k += nv;
3620  }
3621  }
3622  else
3623  {
3624  MPI_Reduce(&NumOfVertices, &TG_nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3625  MPI_Reduce(&NumOfElements, &TG_ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3626  nv = NumOfBdrElements + shared_faces.Size();
3627  MPI_Reduce(&nv, &TG_nbe, 1, MPI_INT, MPI_SUM, 0, MyComm);
3628 
3629  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3630  vert.SetSize(Dim*NumOfVertices);
3631  for (i = 0; i < NumOfVertices; i++)
3632  for (j = 0; j < Dim; j++)
3633  {
3634  vert[Dim*i+j] = vertices[i](j);
3635  }
3636  MPI_Send(&vert[0], Dim*NumOfVertices, MPI_DOUBLE, 0, 445, MyComm);
3637  // elements
3638  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3639  MPI_Send(&NumOfElements, 1, MPI_INT, 0, 446, MyComm);
3640  ints.SetSize(NumOfElements*8);
3641  for (i = 0; i < NumOfElements; i++)
3642  {
3643  v = elements[i]->GetVertices();
3644  for (j = 0; j < 8; j++)
3645  {
3646  ints[8*i+j] = v[j];
3647  }
3648  }
3649  MPI_Send(&ints[0], 8*NumOfElements, MPI_INT, 0, 447, MyComm);
3650  // boundary + shared faces
3651  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3652  ne = NumOfBdrElements + shared_faces.Size();
3653  MPI_Send(&ne, 1, MPI_INT, 0, 446, MyComm);
3654  ints.SetSize(4*ne);
3655  for (i = 0; i < NumOfBdrElements; i++)
3656  {
3657  v = boundary[i]->GetVertices();
3658  for (j = 0; j < 4; j++)
3659  {
3660  ints[4*i+j] = v[j];
3661  }
3662  }
3663  for ( ; i < ne; i++)
3664  {
3665  v = shared_faces[i-NumOfBdrElements]->GetVertices();
3666  for (j = 0; j < 4; j++)
3667  {
3668  ints[4*i+j] = v[j];
3669  }
3670  }
3671  MPI_Send(&ints[0], 4*ne, MPI_INT, 0, 447, MyComm);
3672  }
3673  }
3674 
3675  if (Dim == 2)
3676  {
3677  int i, j, k, attr, nv, ne, p;
3678  Array<int> v;
3679  MPI_Status status;
3680  Array<double> vert;
3681  Array<int> ints;
3682 
3683 
3684  if (MyRank == 0)
3685  {
3686  out << "areamesh2\n\n";
3687 
3688  // print the boundary + shared edges information
3689  nv = NumOfBdrElements + shared_edges.Size();
3690  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3691  out << ne << '\n';
3692  // boundary
3693  for (i = 0; i < NumOfBdrElements; i++)
3694  {
3695  attr = boundary[i]->GetAttribute();
3696  boundary[i]->GetVertices(v);
3697  out << attr << " ";
3698  for (j = 0; j < v.Size(); j++)
3699  {
3700  out << v[j] + 1 << " ";
3701  }
3702  out << '\n';
3703  }
3704  // shared edges
3705  for (i = 0; i < shared_edges.Size(); i++)
3706  {
3707  attr = shared_edges[i]->GetAttribute();
3708  shared_edges[i]->GetVertices(v);
3709  out << attr << " ";
3710  for (j = 0; j < v.Size(); j++)
3711  {
3712  out << v[j] + 1 << " ";
3713  }
3714  out << '\n';
3715  }
3716  k = NumOfVertices;
3717  for (p = 1; p < NRanks; p++)
3718  {
3719  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3720  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3721  ints.SetSize(2*ne);
3722  MPI_Recv(&ints[0], 2*ne, MPI_INT, p, 447, MyComm, &status);
3723  for (i = 0; i < ne; i++)
3724  {
3725  out << p+1;
3726  for (j = 0; j < 2; j++)
3727  {
3728  out << " " << k+ints[i*2+j]+1;
3729  }
3730  out << '\n';
3731  }
3732  k += nv;
3733  }
3734 
3735  // print the elements
3736  nv = NumOfElements;
3737  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3738  out << ne << '\n';
3739  for (i = 0; i < NumOfElements; i++)
3740  {
3741  attr = elements[i]->GetAttribute();
3742  elements[i]->GetVertices(v);
3743  out << 1 << " " << 3 << " ";
3744  for (j = 0; j < v.Size(); j++)
3745  {
3746  out << v[j] + 1 << " ";
3747  }
3748  out << '\n';
3749  }
3750  k = NumOfVertices;
3751  for (p = 1; p < NRanks; p++)
3752  {
3753  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3754  MPI_Recv(&ne, 1, MPI_INT, p, 446, MyComm, &status);
3755  ints.SetSize(3*ne);
3756  MPI_Recv(&ints[0], 3*ne, MPI_INT, p, 447, MyComm, &status);
3757  for (i = 0; i < ne; i++)
3758  {
3759  out << p+1 << " " << 3;
3760  for (j = 0; j < 3; j++)
3761  {
3762  out << " " << k+ints[i*3+j]+1;
3763  }
3764  out << '\n';
3765  }
3766  k += nv;
3767  }
3768 
3769  // print the vertices
3770  ne = NumOfVertices;
3771  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3772  out << nv << '\n';
3773  for (i = 0; i < NumOfVertices; i++)
3774  {
3775  for (j = 0; j < Dim; j++)
3776  {
3777  out << vertices[i](j) << " ";
3778  }
3779  out << '\n';
3780  }
3781  for (p = 1; p < NRanks; p++)
3782  {
3783  MPI_Recv(&nv, 1, MPI_INT, p, 444, MyComm, &status);
3784  vert.SetSize(Dim*nv);
3785  MPI_Recv(&vert[0], Dim*nv, MPI_DOUBLE, p, 445, MyComm, &status);
3786  for (i = 0; i < nv; i++)
3787  {
3788  for (j = 0; j < Dim; j++)
3789  {
3790  out << " " << vert[Dim*i+j];
3791  }
3792  out << '\n';
3793  }
3794  }
3795  }
3796  else
3797  {
3798  // boundary + shared faces
3799  nv = NumOfBdrElements + shared_edges.Size();
3800  MPI_Reduce(&nv, &ne, 1, MPI_INT, MPI_SUM, 0, MyComm);
3801  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3802  ne = NumOfBdrElements + shared_edges.Size();
3803  MPI_Send(&ne, 1, MPI_INT, 0, 446, MyComm);
3804  ints.SetSize(2*ne);
3805  for (i = 0; i < NumOfBdrElements; i++)
3806  {
3807  boundary[i]->GetVertices(v);
3808  for (j = 0; j < 2; j++)
3809  {
3810  ints[2*i+j] = v[j];
3811  }
3812  }
3813  for ( ; i < ne; i++)
3814  {
3815  shared_edges[i-NumOfBdrElements]->GetVertices(v);
3816  for (j = 0; j < 2; j++)
3817  {
3818  ints[2*i+j] = v[j];
3819  }
3820  }
3821  MPI_Send(&ints[0], 2*ne, MPI_INT, 0, 447, MyComm);
3822  // elements
3823  ne = NumOfElements;
3824  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3825  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3826  MPI_Send(&NumOfElements, 1, MPI_INT, 0, 446, MyComm);
3827  ints.SetSize(NumOfElements*3);
3828  for (i = 0; i < NumOfElements; i++)
3829  {
3830  elements[i]->GetVertices(v);
3831  for (j = 0; j < 3; j++)
3832  {
3833  ints[3*i+j] = v[j];
3834  }
3835  }
3836  MPI_Send(&ints[0], 3*NumOfElements, MPI_INT, 0, 447, MyComm);
3837  // vertices
3838  ne = NumOfVertices;
3839  MPI_Reduce(&ne, &nv, 1, MPI_INT, MPI_SUM, 0, MyComm);
3840  MPI_Send(&NumOfVertices, 1, MPI_INT, 0, 444, MyComm);
3841  vert.SetSize(Dim*NumOfVertices);
3842  for (i = 0; i < NumOfVertices; i++)
3843  for (j = 0; j < Dim; j++)
3844  {
3845  vert[Dim*i+j] = vertices[i](j);
3846  }
3847  MPI_Send(&vert[0], Dim*NumOfVertices, MPI_DOUBLE,
3848  0, 445, MyComm);
3849  }
3850  }
3851 }
3852 
3853 void ParMesh::PrintInfo(std::ostream &out)
3854 {
3855  int i;
3856  DenseMatrix J(Dim);
3857  double h_min, h_max, kappa_min, kappa_max, h, kappa;
3858 
3859  if (MyRank == 0)
3860  {
3861  out << "Parallel Mesh Stats:" << '\n';
3862  }
3863 
3864  for (i = 0; i < NumOfElements; i++)
3865  {
3866  GetElementJacobian(i, J);
3867  h = pow(fabs(J.Det()), 1.0/double(Dim));
3868  kappa = J.CalcSingularvalue(0) / J.CalcSingularvalue(Dim-1);
3869  if (i == 0)
3870  {
3871  h_min = h_max = h;
3872  kappa_min = kappa_max = kappa;
3873  }
3874  else
3875  {
3876  if (h < h_min) { h_min = h; }
3877  if (h > h_max) { h_max = h; }
3878  if (kappa < kappa_min) { kappa_min = kappa; }
3879  if (kappa > kappa_max) { kappa_max = kappa; }
3880  }
3881  }
3882 
3883  double gh_min, gh_max, gk_min, gk_max;
3884  MPI_Reduce(&h_min, &gh_min, 1, MPI_DOUBLE, MPI_MIN, 0, MyComm);
3885  MPI_Reduce(&h_max, &gh_max, 1, MPI_DOUBLE, MPI_MAX, 0, MyComm);
3886  MPI_Reduce(&kappa_min, &gk_min, 1, MPI_DOUBLE, MPI_MIN, 0, MyComm);
3887  MPI_Reduce(&kappa_max, &gk_max, 1, MPI_DOUBLE, MPI_MAX, 0, MyComm);
3888 
3889  long ldata[5]; // vert, edge, face, elem, neighbors;
3890  long mindata[5], maxdata[5], sumdata[5];
3891 
3892  // count locally owned vertices, edges, and faces
3893  ldata[0] = GetNV();
3894  ldata[1] = GetNEdges();
3895  ldata[2] = GetNFaces();
3896  ldata[3] = GetNE();
3897  ldata[4] = gtopo.GetNumNeighbors()-1;
3898  for (int gr = 1; gr < GetNGroups(); gr++)
3899  if (!gtopo.IAmMaster(gr)) // we are not the master
3900  {
3901  ldata[0] -= group_svert.RowSize(gr-1);
3902  ldata[1] -= group_sedge.RowSize(gr-1);
3903  ldata[2] -= group_sface.RowSize(gr-1);
3904  }
3905 
3906  MPI_Reduce(ldata, mindata, 5, MPI_LONG, MPI_MIN, 0, MyComm);
3907  MPI_Reduce(ldata, sumdata, 5, MPI_LONG, MPI_SUM, 0, MyComm);
3908  MPI_Reduce(ldata, maxdata, 5, MPI_LONG, MPI_MAX, 0, MyComm);
3909 
3910  if (MyRank == 0)
3911  {
3912  out << '\n'
3913  << " "
3914  << setw(12) << "minimum"
3915  << setw(12) << "average"
3916  << setw(12) << "maximum"
3917  << setw(12) << "total" << '\n';
3918  out << " vertices "
3919  << setw(12) << mindata[0]
3920  << setw(12) << sumdata[0]/NRanks
3921  << setw(12) << maxdata[0]
3922  << setw(12) << sumdata[0] << '\n';
3923  out << " edges "
3924  << setw(12) << mindata[1]
3925  << setw(12) << sumdata[1]/NRanks
3926  << setw(12) << maxdata[1]
3927  << setw(12) << sumdata[1] << '\n';
3928  if (Dim == 3)
3929  out << " faces "
3930  << setw(12) << mindata[2]
3931  << setw(12) << sumdata[2]/NRanks
3932  << setw(12) << maxdata[2]
3933  << setw(12) << sumdata[2] << '\n';
3934  out << " elements "
3935  << setw(12) << mindata[3]
3936  << setw(12) << sumdata[3]/NRanks
3937  << setw(12) << maxdata[3]
3938  << setw(12) << sumdata[3] << '\n';
3939  out << " neighbors "
3940  << setw(12) << mindata[4]
3941  << setw(12) << sumdata[4]/NRanks
3942  << setw(12) << maxdata[4] << '\n';
3943  out << '\n'
3944  << " "
3945  << setw(12) << "minimum"
3946  << setw(12) << "maximum" << '\n';
3947  out << " h "
3948  << setw(12) << gh_min
3949  << setw(12) << gh_max << '\n';
3950  out << " kappa "
3951  << setw(12) << gk_min
3952  << setw(12) << gk_max << '\n';
3953  out << std::flush;
3954  }
3955 }
3956 
3957 long ParMesh::ReduceInt(int value) const
3958 {
3959  long local = value, global;
3960  MPI_Allreduce(&local, &global, 1, MPI_LONG, MPI_SUM, MyComm);
3961  return global;
3962 }
3963 
3965 {
3966  delete pncmesh;
3967  ncmesh = pncmesh = NULL;
3968 
3970 
3971  for (int i = 0; i < shared_faces.Size(); i++)
3972  {
3974  }
3975  for (int i = 0; i < shared_edges.Size(); i++)
3976  {
3978  }
3979 
3980  // The Mesh destructor is called automatically
3981 }
3982 
3983 }
3984 
3985 #endif
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:228
Abstract class for Finite Elements.
Definition: fe.hpp:44
int GetNFaceNeighbors() const
Definition: pmesh.hpp:133
void Create(ListOfIntegerSets &groups, int mpitag)
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:173
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:3218
int Size() const
Logical size of the array.
Definition: array.hpp:109
void Recreate(const int n, const int *p)
Definition: sets.cpp:59
ElementTransformation * Face
Definition: eltrans.hpp:123
int * GetJ()
Definition: table.hpp:108
Class for integration rule.
Definition: intrules.hpp:83
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:529
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
virtual ~ParMesh()
Definition: pmesh.cpp:3964
void SaveAsOne(std::ostream &out=std::cout)
Merge the local grid functions.
Definition: pgridfunc.cpp:327
void UseExternalData(double *ext_data, int i, int j, int k)
Definition: densemat.hpp:606
void FreeElement(Element *E)
Definition: mesh.cpp:8183
void DerefineMesh(const Array< int > &derefinements)
Derefine elements once a list of derefinements is known.
Definition: mesh.cpp:5764
virtual void Update(bool want_transform=true)
Definition: fespace.cpp:1421
int NRanks
Definition: pmesh.hpp:35
virtual Element * Duplicate(Mesh *m) const =0
void GetEdgeVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of edge i.
Definition: mesh.cpp:3434
void AddColumnsInRow(int r, int ncol)
Definition: table.hpp:72
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
Definition: table.cpp:74
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:75
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
void GetFaceElements(int Face, int *Elem1, int *Elem2)
Definition: mesh.cpp:712
Array< Element * > boundary
Definition: mesh.hpp:73
int * GeneratePartitioning(int nparts, int part_method=1)
Definition: mesh.cpp:4338
CoarseFineTransformations CoarseFineTr
Definition: mesh.hpp:115
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:1588
int own_nodes
Definition: mesh.hpp:121
bool Conforming() const
Definition: mesh.hpp:857
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:499
IsoparametricTransformation Transformation
Definition: mesh.hpp:110
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:132
int GetFaceGeometryType(int Face) const
Definition: mesh.cpp:724
Array< Element * > face_nbr_elements
Definition: pmesh.hpp:107
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition: densemat.cpp:416
int NumOfEdges
Definition: mesh.hpp:58
Lists all edges/faces in the nonconforming mesh.
Definition: ncmesh.hpp:167
void GetLocalFaceTransformation(int face_type, int elem_type, IsoparametricTransformation &Transf, int inf)
Used in GetFaceElementTransformations (...)
Definition: mesh.cpp:570
int GetGroupSize(int g) const
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:1547
int BaseGeom
Definition: mesh.hpp:60
Array< int > sface_lface
Definition: pmesh.hpp:48
void SetDims(int rows, int nnz)
Definition: table.cpp:132
const int * GetGroup(int g) const
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:160
T * GetData()
Returns the data.
Definition: array.hpp:91
virtual void LocalRefinement(const Array< int > &marked_el, int type=3)
This function is not public anymore. Use GeneralRefinement instead.
Definition: pmesh.cpp:1632
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
GridFunction * Nodes
Definition: mesh.hpp:120
ElementTransformation * GetGhostFaceTransformation(FaceElementTransformations *FETr, int face_type, int face_geom)
Definition: pmesh.cpp:1432
int NumOfElements
Definition: mesh.hpp:57
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:179
const FiniteElement * GetTraceElement(int i, int geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:1387
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:124
int GetBdrElementEdgeIndex(int i) const
Definition: mesh.cpp:3648
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:496
void SynchronizeDerefinementData(Array< Type > &elem_data, const Table &deref_table)
Definition: pncmesh.cpp:1173
const Element * GetFace(int i) const
Definition: mesh.hpp:547
virtual void SetVertices(const int *ind)
Set the indices the element according to the input.
Definition: element.cpp:17
void GetElementVertices(int i, Array< int > &dofs) const
Returns the indices of the dofs of element i.
Definition: mesh.hpp:558
virtual void OnMeshUpdated(Mesh *mesh)
Definition: pncmesh.cpp:147
Abstract parallel finite element space.
Definition: pfespace.hpp:28
bool have_face_nbr_data
Definition: pmesh.hpp:103
const FiniteElement * GetFE() const
Definition: eltrans.hpp:82
Array< int > face_nbr_vertices_offset
Definition: pmesh.hpp:106
std::vector< MeshId > conforming
Definition: ncmesh.hpp:169
Data type for vertex.
Definition: vertex.hpp:21
Array< int > face_nbr_group
Definition: pmesh.hpp:104
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3212
Data type quadrilateral element.
void GetVertices(Vector &vert_coord) const
Definition: mesh.cpp:5043
int BaseBdrGeom
Definition: mesh.hpp:60
ElementTransformation * Elem2
Definition: eltrans.hpp:123
int GetFaceElementType(int Face) const
Definition: mesh.cpp:729
void CheckBdrElementOrientation(bool fix_it=true)
Check the orientation of the boundary elements.
Definition: mesh.cpp:3260
int Size_of_connections() const
Definition: table.hpp:92
Array< Element * > faces
Definition: mesh.hpp:74
const IntegrationRule * GetVertices(int GeomType)
Definition: geom.cpp:172
const NCList & GetSharedEdges()
Definition: pncmesh.hpp:106
void GetVertexToVertexTable(DSTable &) const
Definition: mesh.cpp:3759
Array< int > sedge_ledge
Definition: pmesh.hpp:47
void RedRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:215
Operation last_operation
Definition: mesh.hpp:147
int GetGeometryType() const
Definition: element.hpp:47
void DeleteAll()
Delete whole array.
Definition: array.hpp:463
void AddConnections(int r, const int *c, int nc)
Definition: table.cpp:96
void InitRefinementTransforms()
Definition: mesh.cpp:6492
IsoparametricTransformation FaceTransformation
Definition: mesh.hpp:111
Array< NCFaceInfo > nc_faces_info
Definition: mesh.hpp:99
bool IAmMaster(int g) const
bool WantSkipSharedMaster(const NCMesh::Master &master) const
Definition: pmesh.cpp:2914
A parallel extension of the NCMesh class.
Definition: pncmesh.hpp:65
Element * NewElement(int geom)
Definition: mesh.cpp:2193
Table * el_to_face
Definition: mesh.hpp:102
FaceElementTransformations FaceElemTr
Definition: mesh.hpp:112
ParNCMesh * pncmesh
Definition: pmesh.hpp:113
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:3039
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:801
void ExchangeFaceNbrData()
Definition: pmesh.cpp:896
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2116
void PrintAsOne(std::ostream &out=std::cout)
Definition: pmesh.cpp:3055
virtual void GetVertices(Array< int > &v) const
Returns the indices of the element&#39;s vertices.
Definition: segment.cpp:40
void GetElementJacobian(int i, DenseMatrix &J)
Definition: mesh.cpp:35
Geometry Geometries
Definition: geom.cpp:544
void Rebalance()
Load balance the mesh. NC meshes only.
Definition: pmesh.cpp:2255
bool Nonconforming() const
Definition: mesh.hpp:858
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:835
Array< int > face_nbr_elements_offset
Definition: pmesh.hpp:105
MPI_Comm MyComm
Definition: pmesh.hpp:34
Symmetric 3D Table.
Definition: stable3d.hpp:28
const Table & GetDerefinementTable()
Definition: ncmesh.cpp:1356
DenseTensor point_matrices
matrices for IsoparametricTransformation
Definition: ncmesh.hpp:53
Array< Element * > shared_edges
Definition: pmesh.hpp:37
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:376
virtual void NonconformingRefinement(const Array< Refinement > &refinements, int nc_limit=0)
This function is not public anymore. Use GeneralRefinement instead.
Definition: pmesh.cpp:2155
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:124
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition: pmesh.cpp:1566
int GetNumNeighbors() const
void Clear()
Definition: table.cpp:346
void GetFaceTransformation(int i, IsoparametricTransformation *FTr)
Definition: mesh.cpp:344
void CreateRefinementFlag(int refinement_edges[2], int type, int flag=0)
Definition: tetrahedron.cpp:59
std::vector< Master > masters
Definition: ncmesh.hpp:170
static FiniteElement * GetTransformationFEforElementType(int)
Definition: mesh.cpp:221
void AddConnection(int r, int c)
Definition: table.hpp:74
STable3D * GetElementToFaceTable(int ret_ftbl=0)
Definition: mesh.cpp:4111
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:1462
static void Rotate3(int &, int &, int &)
Definition: mesh.hpp:991
int MeshGenerator()
Definition: mesh.hpp:489
Table send_face_nbr_vertices
Definition: pmesh.hpp:111
int GetFaceSplittings(Element *face, const DSTable &v_to_v, int *middle)
Return a number(0-4) identifying how the given face has been split.
Definition: pmesh.cpp:744
T Max() const
Definition: array.cpp:90
void InitTables()
Definition: mesh.cpp:745
int InitialPartition(int index) const
Helper to get the partition when the serial mesh is being split initially.
Definition: pncmesh.hpp:278
void CheckElementOrientation(bool fix_it=true)
Check the orientation of the elements.
Definition: mesh.cpp:3054
const IntegrationRule & GetNodes() const
Definition: fe.hpp:113
int Insert(IntegerSet &s)
Definition: sets.cpp:82
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition: array.hpp:122
Array< Embedding > embeddings
fine element positions in their parents
Definition: ncmesh.hpp:54
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:391
Array< Element * > shared_faces
Definition: pmesh.hpp:38
int GetAttribute() const
Return element&#39;s attribute.
Definition: element.hpp:50
int GetElementToEdgeTable(Table &, Array< int > &)
Definition: mesh.cpp:3784
Data type triangle element.
Definition: triangle.hpp:23
int GetElementType(int i) const
Returns the type of element i.
Definition: mesh.cpp:3680
void MarkCoarseLevel()
Definition: ncmesh.cpp:2732
const Element * GetElement(int i) const
Definition: mesh.hpp:539
virtual void QuadUniformRefinement()
Refine quadrilateral mesh.
Definition: pmesh.cpp:2475
IsoparametricTransformation Transformation2
Definition: mesh.hpp:110
void Init()
Definition: mesh.cpp:734
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:86
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:151
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:263
virtual void ReorientTetMesh()
Definition: mesh.cpp:4184
virtual void HexUniformRefinement()
Refine a hexahedral mesh.
Definition: pmesh.cpp:2562
int NumOfBdrElements
Definition: mesh.hpp:57
Table * el_to_edge
Definition: mesh.hpp:101
void Rebalance()
Definition: pncmesh.cpp:1309
int slaves_end
slave faces
Definition: ncmesh.hpp:145
Data type tetrahedron element.
Definition: tetrahedron.hpp:22
Vector & FaceNbrData()
Definition: pgridfunc.hpp:115
void GetBdrElementFace(int i, int *, int *) const
Return the index and the orientation of the face of bdr element i. (3D)
Definition: mesh.cpp:3603
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
Definition: mesh.cpp:6736
Array< int > bdr_attributes
Definition: mesh.hpp:141
void RefineGroups(const DSTable &v_to_v, int *middle)
Update the groups after tet refinement.
Definition: pmesh.cpp:2288
virtual bool NonconformingDerefinement(Array< double > &elem_error, double threshold, int nc_limit=0, int op=1)
NC version of GeneralDerefinement.
Definition: pmesh.cpp:2207
Array< Vertex > face_nbr_vertices
Definition: pmesh.hpp:108
std::vector< Slave > slaves
Definition: ncmesh.hpp:171
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
void ApplyLocalSlaveTransformation(IsoparametricTransformation &transf, const FaceInfo &fi)
Definition: mesh.cpp:673
static void PrintElement(const Element *, std::ostream &)
Definition: mesh.cpp:2254
int GetNeighborRank(int i) const
void AddAColumnInRow(int r)
Definition: table.hpp:71
void mfem_error(const char *msg)
Definition: error.cpp:23
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:331
void SetSubVector(const Array< int > &dofs, const Vector &elemvect)
Definition: vector.cpp:493
virtual long ReduceInt(int value) const
Utility function: sum integers from all processors (Allreduce).
Definition: pmesh.cpp:3957
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:243
bool IsSlaveFace(const FaceInfo &fi) const
Definition: mesh.cpp:668
Array< Vertex > vertices
Definition: mesh.hpp:72
virtual void PrintInfo(std::ostream &out=std::cout)
Print various parallel mesh stats.
Definition: pmesh.cpp:3853
virtual void QuadUniformRefinement()
Refine quadrilateral mesh.
Definition: mesh.cpp:5192
virtual void HexUniformRefinement()
Refine hexahedral mesh.
Definition: mesh.cpp:5294
void Swap(Mesh &other, bool non_geometry=false)
Definition: mesh.cpp:5922
Array< Element * > elements
Definition: mesh.hpp:69
void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:141
void ShiftUpI()
Definition: table.cpp:107
int meshgen
Definition: mesh.hpp:62
virtual void PrintXG(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2719
int GetElementBaseGeometry(int i=0) const
Definition: mesh.hpp:551
static FiniteElementCollection * New(const char *name)
Definition: fe_coll.cpp:44
virtual void GetVertices(Array< int > &v) const
Returns the indices of the element&#39;s vertices.
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:824
NURBSExtension * NURBSext
Definition: mesh.hpp:143
int GetNV() const
Definition: mesh.hpp:493
Array< int > be_to_edge
Definition: mesh.hpp:104
virtual const char * Name() const
Definition: fe_coll.hpp:50
virtual void Refine(const Array< Refinement > &refinements)
Definition: pncmesh.cpp:911
void GroupFace(int group, int i, int &face, int &o)
Definition: pmesh.cpp:705
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:237
Table group_sedge
Definition: pmesh.hpp:42
void ExchangeFaceNbrNodes()
Definition: pmesh.cpp:1299
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:144
A set of integers.
Definition: sets.hpp:23
void SetIJ(int *newI, int *newJ, int newsize=-1)
Definition: table.cpp:196
void MakeJ()
Definition: table.cpp:84
Table group_svert
Shared objects in each group.
Definition: pmesh.hpp:41
int GroupNFaces(int group)
Definition: pmesh.hpp:120
void SetFE(const FiniteElement *FE)
Definition: eltrans.hpp:81
long sequence
Definition: mesh.hpp:67
IsoparametricTransformation Transf
Definition: eltrans.hpp:114
ElementTransformation * Elem1
Definition: eltrans.hpp:123
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1652
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:429
Array< FaceInfo > faces_info
Definition: mesh.hpp:98
virtual int GetNVertices() const =0
void SetAttribute(const int attr)
Set element&#39;s attribute.
Definition: element.hpp:53
STable3D * GetFacesTable()
Definition: mesh.cpp:4076
NCMesh * ncmesh
Definition: mesh.hpp:144
int Dim
Definition: mesh.hpp:54
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
Definition: pncmesh.cpp:1253
int GetEdgeSplittings(Element *edge, const DSTable &v_to_v, int *middle)
Return a number(0-1) identifying how the given edge has been split.
Definition: pmesh.cpp:722
int GetNEdges() const
Return the number of edges.
Definition: mesh.hpp:502
void GroupEdge(int group, int i, int &edge, int &o)
Definition: pmesh.cpp:697
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:505
const NCList & GetSharedList(int type)
Helper to get shared vertices/edges/faces (&#39;type&#39; == 0/1/2 resp.).
Definition: pncmesh.hpp:121
double kappa
Definition: ex3.cpp:46
Vector data type.
Definition: vector.hpp:33
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:175
virtual void Transform(const IntegrationPoint &, Vector &)=0
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5114
void AverageVertices(int *indexes, int n, int result)
Definition: mesh.cpp:5162
void GenerateOffsets(int N, HYPRE_Int loc_sizes[], Array< HYPRE_Int > *offsets[]) const
Definition: pmesh.cpp:781
int NumberOfEntries() const
Definition: table.hpp:226
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of &quot;test&quot; relative to &quot;base&quot;.
Definition: mesh.cpp:3164
bool IsGhost(int type, int index) const
Returns true if the specified vertex/edge/face is a ghost.
Definition: pncmesh.hpp:177
int * GetI()
Definition: table.hpp:107
void PrintAsOneXG(std::ostream &out=std::cout)
Old mesh format (Netgen/Truegrid) version of &#39;PrintAsOne&#39;.
Definition: pmesh.cpp:3330
Array< int > svert_lvert
Shared to local index mapping.
Definition: pmesh.hpp:46
void GenerateFaces()
Definition: mesh.cpp:3960
const NCList & GetSharedFaces()
Definition: pncmesh.hpp:114
int MyRank
Definition: pmesh.hpp:35
int spaceDim
Definition: mesh.hpp:55
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Definition: ncmesh.hpp:188
void ParseRefinementFlag(int refinement_edges[2], int &type, int &flag)
Definition: tetrahedron.cpp:43
int RowSize(int i) const
Definition: table.hpp:102
Class for parallel grid function.
Definition: pgridfunc.hpp:31
friend class ParNCMesh
Definition: mesh.hpp:49
Table group_sface
Definition: pmesh.hpp:43
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:71
Table send_face_nbr_elements
Definition: pmesh.hpp:110
List of integer sets.
Definition: sets.hpp:49
int NumOfVertices
Definition: mesh.hpp:57
Table * GetFaceToAllElementTable() const
Definition: pmesh.cpp:1374
void DeleteFaceNbrData()
Definition: pmesh.cpp:875
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
Definition: mesh.cpp:5864
GroupTopology gtopo
Definition: pmesh.hpp:100
void Bisection(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6159
int GetNGroups()
Definition: pmesh.hpp:115
Class for parallel meshes.
Definition: pmesh.hpp:28
Abstract data type element.
Definition: element.hpp:27
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:716
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:1357
const Table & ElementToEdgeTable() const
Definition: mesh.cpp:3862
Data type line segment element.
Definition: segment.hpp:22
void GenerateNCFaceInfo()
Definition: mesh.cpp:4029
virtual void LimitNCLevel(int max_nc_level)
Parallel version of NCMesh::LimitNCLevel.
Definition: pncmesh.cpp:983
Array< int > attributes
Definition: mesh.hpp:140
Table * GetVertexToElementTable()
The returned Table must be destroyed by the caller.
Definition: mesh.cpp:3496
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:543
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:123
virtual int GetType() const =0
Returns element&#39;s type.
void UpdateNodes()
Update the nodes of a curved mesh after refinement.
Definition: mesh.cpp:5183
Defines the position of a fine element within a coarse element.
Definition: ncmesh.hpp:42
void GetFaceNeighbors(ParMesh &pmesh)
Definition: pncmesh.cpp:594
int GroupNEdges(int group)
Definition: pmesh.hpp:119
void GreenRefinement(int i, const DSTable &v_to_v, int *edge1, int *edge2, int *middle)
Definition: mesh.hpp:221
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:61
virtual void NURBSUniformRefinement()
Refine NURBS mesh.
Definition: pmesh.cpp:2711
virtual void Print(std::ostream &out=std::cout) const
Definition: pmesh.cpp:2927
int NumOfFaces
Definition: mesh.hpp:58