MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pfespace.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "pfespace.hpp"
17 #include "prestriction.hpp"
18 #include "../general/forall.hpp"
19 #include "../general/sort_pairs.hpp"
20 #include "../mesh/mesh_headers.hpp"
21 #include "../general/binaryio.hpp"
22 
23 #include <climits> // INT_MAX
24 #include <limits>
25 #include <list>
26 
27 namespace mfem
28 {
29 
31  const ParFiniteElementSpace &orig, ParMesh *pmesh,
32  const FiniteElementCollection *fec)
33  : FiniteElementSpace(orig, pmesh, fec)
34 {
35  ParInit(pmesh ? pmesh : orig.pmesh);
36 }
37 
39  const FiniteElementSpace &orig, ParMesh &pmesh,
40  const FiniteElementCollection *fec)
41  : FiniteElementSpace(orig, &pmesh, fec)
42 {
43  ParInit(&pmesh);
44 }
45 
47  ParMesh *pm, const FiniteElementSpace *global_fes, const int *partitioning,
49  : FiniteElementSpace(pm, MakeLocalNURBSext(global_fes->GetNURBSext(),
50  pm->NURBSext),
51  f ? f : global_fes->FEColl(),
52  global_fes->GetVDim(), global_fes->GetOrdering())
53 {
54  ParInit(pm);
55  // For NURBS spaces, the variable-order data is contained in the
56  // NURBSExtension of 'global_fes' and inside the ParNURBSExtension of 'pm'.
57 
58  // TODO: when general variable-order support is added, copy the local portion
59  // of the variable-order data from 'global_fes' to 'this'.
60 }
61 
63  ParMesh *pm, const FiniteElementCollection *f, int dim, int ordering)
64  : FiniteElementSpace(pm, f, dim, ordering)
65 {
66  ParInit(pm);
67 }
68 
71  int dim, int ordering)
72  : FiniteElementSpace(pm, ext, f, dim, ordering)
73 {
74  ParInit(pm);
75 }
76 
77 // static method
78 ParNURBSExtension *ParFiniteElementSpace::MakeLocalNURBSext(
79  const NURBSExtension *globNURBSext, const NURBSExtension *parNURBSext)
80 {
81  if (globNURBSext == NULL) { return NULL; }
82  const ParNURBSExtension *pNURBSext =
83  dynamic_cast<const ParNURBSExtension*>(parNURBSext);
84  MFEM_ASSERT(pNURBSext, "need a ParNURBSExtension");
85  // make a copy of globNURBSext:
86  NURBSExtension *tmp_globNURBSext = new NURBSExtension(*globNURBSext);
87  // tmp_globNURBSext will be deleted by the following ParNURBSExtension ctor:
88  return new ParNURBSExtension(tmp_globNURBSext, pNURBSext);
89 }
90 
91 void ParFiniteElementSpace::ParInit(ParMesh *pm)
92 {
93  pmesh = pm;
94  pncmesh = pm->pncmesh;
95 
96  MyComm = pmesh->GetComm();
97  NRanks = pmesh->GetNRanks();
98  MyRank = pmesh->GetMyRank();
99 
100  gcomm = NULL;
101 
102  P = NULL;
103  Pconf = NULL;
104  R = NULL;
105 
106  num_face_nbr_dofs = -1;
107 
108  if (NURBSext && !pNURBSext())
109  {
110  // This is necessary in some cases: e.g. when the FiniteElementSpace
111  // constructor creates a serial NURBSExtension of higher order than the
112  // mesh NURBSExtension.
113  MFEM_ASSERT(own_ext, "internal error");
114 
115  ParNURBSExtension *pNe = new ParNURBSExtension(
116  NURBSext, dynamic_cast<ParNURBSExtension *>(pmesh->NURBSext));
117  // serial NURBSext is destroyed by the above constructor
118  NURBSext = pNe;
119  UpdateNURBS();
120  }
121 
122  Construct(); // parallel version of Construct().
123 
124  // Apply the ldof_signs to the elem_dof Table
125  if (Conforming() && !NURBSext)
126  {
127  ApplyLDofSigns(*elem_dof);
128  }
129 }
130 
131 void ParFiniteElementSpace::Construct()
132 {
133  if (NURBSext)
134  {
135  ConstructTrueNURBSDofs();
136  GenerateGlobalOffsets();
137  }
138  else if (Conforming())
139  {
140  ConstructTrueDofs();
141  GenerateGlobalOffsets();
142  }
143  else // Nonconforming()
144  {
145  // Initialize 'gcomm' for the cut (aka "partially conforming") space.
146  // In the process, the array 'ldof_ltdof' is also initialized (for the cut
147  // space) and used; however, it will be overwritten below with the real
148  // true dofs. Also, 'ldof_sign' and 'ldof_group' are constructed for the
149  // cut space.
150  ConstructTrueDofs();
151 
152  ngedofs = ngfdofs = 0;
153 
154  // calculate number of ghost DOFs
155  ngvdofs = pncmesh->GetNGhostVertices()
157 
158  if (pmesh->Dimension() > 1)
159  {
160  ngedofs = pncmesh->GetNGhostEdges()
162  }
163 
164  if (pmesh->Dimension() > 2)
165  {
166  int stride = fec->DofForGeometry(Geometry::SQUARE);
167  ngfdofs = pncmesh->GetNGhostFaces() * stride;
168  }
169 
170  // total number of ghost DOFs. Ghost DOFs start at index 'ndofs', i.e.,
171  // after all regular DOFs
172  ngdofs = ngvdofs + ngedofs + ngfdofs;
173 
174  // get P and R matrices, initialize DOF offsets, etc. NOTE: in the NC
175  // case this needs to be done here to get the number of true DOFs
176  ltdof_size = BuildParallelConformingInterpolation(
177  &P, &R, dof_offsets, tdof_offsets, &ldof_ltdof, false);
178 
179  // TODO future: split BuildParallelConformingInterpolation into two parts
180  // to overlap its communication with processing between this constructor
181  // and the point where the P matrix is actually needed.
182  }
183 }
184 
186 {
187  long ltdofs = ltdof_size;
188  long min_ltdofs, max_ltdofs, sum_ltdofs;
189 
190  MPI_Reduce(&ltdofs, &min_ltdofs, 1, MPI_LONG, MPI_MIN, 0, MyComm);
191  MPI_Reduce(&ltdofs, &max_ltdofs, 1, MPI_LONG, MPI_MAX, 0, MyComm);
192  MPI_Reduce(&ltdofs, &sum_ltdofs, 1, MPI_LONG, MPI_SUM, 0, MyComm);
193 
194  if (MyRank == 0)
195  {
196  double avg = double(sum_ltdofs) / NRanks;
197  mfem::out << "True DOF partitioning: min " << min_ltdofs
198  << ", avg " << std::fixed << std::setprecision(1) << avg
199  << ", max " << max_ltdofs
200  << ", (max-avg)/avg " << 100.0*(max_ltdofs - avg)/avg
201  << "%" << std::endl;
202  }
203 
204  if (NRanks <= 32)
205  {
206  if (MyRank == 0)
207  {
208  mfem::out << "True DOFs by rank: " << ltdofs;
209  for (int i = 1; i < NRanks; i++)
210  {
211  MPI_Status status;
212  MPI_Recv(&ltdofs, 1, MPI_LONG, i, 123, MyComm, &status);
213  mfem::out << " " << ltdofs;
214  }
215  mfem::out << "\n";
216  }
217  else
218  {
219  MPI_Send(&ltdofs, 1, MPI_LONG, 0, 123, MyComm);
220  }
221  }
222 }
223 
224 void ParFiniteElementSpace::GetGroupComm(
225  GroupCommunicator &gc, int ldof_type, Array<int> *ldof_sign)
226 {
227  int gr;
228  int ng = pmesh->GetNGroups();
229  int nvd, ned, ntd = 0, nqd = 0;
230  Array<int> dofs;
231 
232  int group_ldof_counter;
233  Table &group_ldof = gc.GroupLDofTable();
234 
237 
238  if (fdofs)
239  {
241  {
243  }
245  {
247  }
248  }
249 
250  if (ldof_sign)
251  {
252  ldof_sign->SetSize(GetNDofs());
253  *ldof_sign = 1;
254  }
255 
256  // count the number of ldofs in all groups (excluding the local group 0)
257  group_ldof_counter = 0;
258  for (gr = 1; gr < ng; gr++)
259  {
260  group_ldof_counter += nvd * pmesh->GroupNVertices(gr);
261  group_ldof_counter += ned * pmesh->GroupNEdges(gr);
262  group_ldof_counter += ntd * pmesh->GroupNTriangles(gr);
263  group_ldof_counter += nqd * pmesh->GroupNQuadrilaterals(gr);
264  }
265  if (ldof_type)
266  {
267  group_ldof_counter *= vdim;
268  }
269  // allocate the I and J arrays in group_ldof
270  group_ldof.SetDims(ng, group_ldof_counter);
271 
272  // build the full group_ldof table
273  group_ldof_counter = 0;
274  group_ldof.GetI()[0] = group_ldof.GetI()[1] = 0;
275  for (gr = 1; gr < ng; gr++)
276  {
277  int j, k, l, m, o, nv, ne, nt, nq;
278  const int *ind;
279 
280  nv = pmesh->GroupNVertices(gr);
281  ne = pmesh->GroupNEdges(gr);
282  nt = pmesh->GroupNTriangles(gr);
283  nq = pmesh->GroupNQuadrilaterals(gr);
284 
285  // vertices
286  if (nvd > 0)
287  {
288  for (j = 0; j < nv; j++)
289  {
290  k = pmesh->GroupVertex(gr, j);
291 
292  dofs.SetSize(nvd);
293  m = nvd * k;
294  for (l = 0; l < nvd; l++, m++)
295  {
296  dofs[l] = m;
297  }
298 
299  if (ldof_type)
300  {
301  DofsToVDofs(dofs);
302  }
303 
304  for (l = 0; l < dofs.Size(); l++)
305  {
306  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
307  }
308  }
309  }
310 
311  // edges
312  if (ned > 0)
313  {
314  for (j = 0; j < ne; j++)
315  {
316  pmesh->GroupEdge(gr, j, k, o);
317 
318  dofs.SetSize(ned);
319  m = nvdofs+k*ned;
321  for (l = 0; l < ned; l++)
322  {
323  if (ind[l] < 0)
324  {
325  dofs[l] = m + (-1-ind[l]);
326  if (ldof_sign)
327  {
328  (*ldof_sign)[dofs[l]] = -1;
329  }
330  }
331  else
332  {
333  dofs[l] = m + ind[l];
334  }
335  }
336 
337  if (ldof_type)
338  {
339  DofsToVDofs(dofs);
340  }
341 
342  for (l = 0; l < dofs.Size(); l++)
343  {
344  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
345  }
346  }
347  }
348 
349  // triangles
350  if (ntd > 0)
351  {
352  for (j = 0; j < nt; j++)
353  {
354  pmesh->GroupTriangle(gr, j, k, o);
355 
356  dofs.SetSize(ntd);
357  m = nvdofs+nedofs+fdofs[k];
359  for (l = 0; l < ntd; l++)
360  {
361  if (ind[l] < 0)
362  {
363  dofs[l] = m + (-1-ind[l]);
364  if (ldof_sign)
365  {
366  (*ldof_sign)[dofs[l]] = -1;
367  }
368  }
369  else
370  {
371  dofs[l] = m + ind[l];
372  }
373  }
374 
375  if (ldof_type)
376  {
377  DofsToVDofs(dofs);
378  }
379 
380  for (l = 0; l < dofs.Size(); l++)
381  {
382  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
383  }
384  }
385  }
386 
387  // quadrilaterals
388  if (nqd > 0)
389  {
390  for (j = 0; j < nq; j++)
391  {
392  pmesh->GroupQuadrilateral(gr, j, k, o);
393 
394  dofs.SetSize(nqd);
395  m = nvdofs+nedofs+fdofs[k];
397  for (l = 0; l < nqd; l++)
398  {
399  if (ind[l] < 0)
400  {
401  dofs[l] = m + (-1-ind[l]);
402  if (ldof_sign)
403  {
404  (*ldof_sign)[dofs[l]] = -1;
405  }
406  }
407  else
408  {
409  dofs[l] = m + ind[l];
410  }
411  }
412 
413  if (ldof_type)
414  {
415  DofsToVDofs(dofs);
416  }
417 
418  for (l = 0; l < dofs.Size(); l++)
419  {
420  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
421  }
422  }
423  }
424 
425  group_ldof.GetI()[gr+1] = group_ldof_counter;
426  }
427 
428  gc.Finalize();
429 }
430 
431 void ParFiniteElementSpace::ApplyLDofSigns(Array<int> &dofs) const
432 {
433  MFEM_ASSERT(Conforming(), "wrong code path");
434 
435  for (int i = 0; i < dofs.Size(); i++)
436  {
437  if (dofs[i] < 0)
438  {
439  if (ldof_sign[-1-dofs[i]] < 0)
440  {
441  dofs[i] = -1-dofs[i];
442  }
443  }
444  else
445  {
446  if (ldof_sign[dofs[i]] < 0)
447  {
448  dofs[i] = -1-dofs[i];
449  }
450  }
451  }
452 }
453 
454 void ParFiniteElementSpace::ApplyLDofSigns(Table &el_dof) const
455 {
456  Array<int> all_dofs(el_dof.GetJ(), el_dof.Size_of_connections());
457  ApplyLDofSigns(all_dofs);
458 }
459 
461 {
462  if (elem_dof)
463  {
464  elem_dof->GetRow(i, dofs);
465  return;
466  }
468  if (Conforming())
469  {
470  ApplyLDofSigns(dofs);
471  }
472 }
473 
475 {
476  if (bdrElem_dof)
477  {
478  bdrElem_dof->GetRow(i, dofs);
479  return;
480  }
482  if (Conforming())
483  {
484  ApplyLDofSigns(dofs);
485  }
486 }
487 
489 {
490  if (face_dof)
491  {
492  face_dof->GetRow(i, dofs);
493  return;
494  }
496  if (Conforming())
497  {
498  ApplyLDofSigns(dofs);
499  }
500 }
501 
502 
504  ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul) const
505 {
506  const bool is_dg_space = IsDGSpace();
507  const L2FaceValues m = (is_dg_space && mul==L2FaceValues::DoubleValued) ?
509  auto key = std::make_tuple(is_dg_space, e_ordering, type, m);
510  auto itr = L2F.find(key);
511  if (itr != L2F.end())
512  {
513  return itr->second;
514  }
515  else
516  {
517  Operator* res;
518  if (is_dg_space)
519  {
520  res = new ParL2FaceRestriction(*this, e_ordering, type, m);
521  }
522  else
523  {
524  res = new H1FaceRestriction(*this, e_ordering, type);
525  }
526  L2F[key] = res;
527  return res;
528  }
529 }
530 
532  int group, int ei, Array<int> &dofs) const
533 {
534  int l_edge, ori;
535  MFEM_ASSERT(0 <= ei && ei < pmesh->GroupNEdges(group), "invalid edge index");
536  pmesh->GroupEdge(group, ei, l_edge, ori);
537  if (ori > 0) // ori = +1 or -1
538  {
539  GetEdgeDofs(l_edge, dofs);
540  }
541  else
542  {
543  Array<int> rdofs;
544  fec->SubDofOrder(Geometry::SEGMENT, 1, 1, dofs);
545  GetEdgeDofs(l_edge, rdofs);
546  for (int i = 0; i < dofs.Size(); i++)
547  {
548  const int di = dofs[i];
549  dofs[i] = (di >= 0) ? rdofs[di] : -1-rdofs[-1-di];
550  }
551  }
552 }
553 
555  int group, int fi, Array<int> &dofs) const
556 {
557  int l_face, ori;
558  MFEM_ASSERT(0 <= fi && fi < pmesh->GroupNTriangles(group),
559  "invalid triangular face index");
560  pmesh->GroupTriangle(group, fi, l_face, ori);
561  if (ori == 0)
562  {
563  GetFaceDofs(l_face, dofs);
564  }
565  else
566  {
567  Array<int> rdofs;
568  fec->SubDofOrder(Geometry::TRIANGLE, 2, ori, dofs);
569  GetFaceDofs(l_face, rdofs);
570  for (int i = 0; i < dofs.Size(); i++)
571  {
572  const int di = dofs[i];
573  dofs[i] = (di >= 0) ? rdofs[di] : -1-rdofs[-1-di];
574  }
575  }
576 }
577 
579  int group, int fi, Array<int> &dofs) const
580 {
581  int l_face, ori;
582  MFEM_ASSERT(0 <= fi && fi < pmesh->GroupNQuadrilaterals(group),
583  "invalid quadrilateral face index");
584  pmesh->GroupQuadrilateral(group, fi, l_face, ori);
585  if (ori == 0)
586  {
587  GetFaceDofs(l_face, dofs);
588  }
589  else
590  {
591  Array<int> rdofs;
592  fec->SubDofOrder(Geometry::SQUARE, 2, ori, dofs);
593  GetFaceDofs(l_face, rdofs);
594  for (int i = 0; i < dofs.Size(); i++)
595  {
596  const int di = dofs[i];
597  dofs[i] = (di >= 0) ? rdofs[di] : -1-rdofs[-1-di];
598  }
599  }
600 }
601 
602 void ParFiniteElementSpace::GenerateGlobalOffsets() const
603 {
604  MFEM_ASSERT(Conforming(), "wrong code path");
605 
606  HYPRE_Int ldof[2];
607  Array<HYPRE_Int> *offsets[2] = { &dof_offsets, &tdof_offsets };
608 
609  ldof[0] = GetVSize();
610  ldof[1] = TrueVSize();
611 
612  pmesh->GenerateOffsets(2, ldof, offsets);
613 
614  if (HYPRE_AssumedPartitionCheck())
615  {
616  // communicate the neighbor offsets in tdof_nb_offsets
617  GroupTopology &gt = GetGroupTopo();
618  int nsize = gt.GetNumNeighbors()-1;
619  MPI_Request *requests = new MPI_Request[2*nsize];
620  MPI_Status *statuses = new MPI_Status[2*nsize];
621  tdof_nb_offsets.SetSize(nsize+1);
622  tdof_nb_offsets[0] = tdof_offsets[0];
623 
624  // send and receive neighbors' local tdof offsets
625  int request_counter = 0;
626  for (int i = 1; i <= nsize; i++)
627  {
628  MPI_Irecv(&tdof_nb_offsets[i], 1, HYPRE_MPI_INT,
629  gt.GetNeighborRank(i), 5365, MyComm,
630  &requests[request_counter++]);
631  }
632  for (int i = 1; i <= nsize; i++)
633  {
634  MPI_Isend(&tdof_nb_offsets[0], 1, HYPRE_MPI_INT,
635  gt.GetNeighborRank(i), 5365, MyComm,
636  &requests[request_counter++]);
637  }
638  MPI_Waitall(request_counter, requests, statuses);
639 
640  delete [] statuses;
641  delete [] requests;
642  }
643 }
644 
645 void ParFiniteElementSpace::Build_Dof_TrueDof_Matrix() const // matrix P
646 {
647  MFEM_ASSERT(Conforming(), "wrong code path");
648 
649  if (P) { return; }
650 
651  int ldof = GetVSize();
652  int ltdof = TrueVSize();
653 
654  HYPRE_Int *i_diag = Memory<HYPRE_Int>(ldof+1);
655  HYPRE_Int *j_diag = Memory<HYPRE_Int>(ltdof);
656  int diag_counter;
657 
658  HYPRE_Int *i_offd = Memory<HYPRE_Int>(ldof+1);
659  HYPRE_Int *j_offd = Memory<HYPRE_Int>(ldof-ltdof);
660  int offd_counter;
661 
662  HYPRE_Int *cmap = Memory<HYPRE_Int>(ldof-ltdof);
663 
664  HYPRE_Int *col_starts = GetTrueDofOffsets();
665  HYPRE_Int *row_starts = GetDofOffsets();
666 
667  Array<Pair<HYPRE_Int, int> > cmap_j_offd(ldof-ltdof);
668 
669  i_diag[0] = i_offd[0] = 0;
670  diag_counter = offd_counter = 0;
671  for (int i = 0; i < ldof; i++)
672  {
673  int ltdof = GetLocalTDofNumber(i);
674  if (ltdof >= 0)
675  {
676  j_diag[diag_counter++] = ltdof;
677  }
678  else
679  {
680  cmap_j_offd[offd_counter].one = GetGlobalTDofNumber(i);
681  cmap_j_offd[offd_counter].two = offd_counter;
682  offd_counter++;
683  }
684  i_diag[i+1] = diag_counter;
685  i_offd[i+1] = offd_counter;
686  }
687 
688  SortPairs<HYPRE_Int, int>(cmap_j_offd, offd_counter);
689 
690  for (int i = 0; i < offd_counter; i++)
691  {
692  cmap[i] = cmap_j_offd[i].one;
693  j_offd[cmap_j_offd[i].two] = i;
694  }
695 
696  P = new HypreParMatrix(MyComm, MyRank, NRanks, row_starts, col_starts,
697  i_diag, j_diag, i_offd, j_offd, cmap, offd_counter);
698 
699  SparseMatrix Pdiag;
700  P->GetDiag(Pdiag);
701  R = Transpose(Pdiag);
702 }
703 
705 {
706  HypreParMatrix *P_pc;
707  Array<HYPRE_Int> P_pc_row_starts, P_pc_col_starts;
708  BuildParallelConformingInterpolation(&P_pc, NULL, P_pc_row_starts,
709  P_pc_col_starts, NULL, true);
710  P_pc->CopyRowStarts();
711  P_pc->CopyColStarts();
712  return P_pc;
713 }
714 
716 {
717  GroupTopology &gt = GetGroupTopo();
718  for (int i = 0; i < ldof_group.Size(); i++)
719  {
720  if (gt.IAmMaster(ldof_group[i])) // we are the master
721  {
722  if (ldof_ltdof[i] >= 0) // see note below
723  {
724  vec[ldof_ltdof[i]] /= gt.GetGroupSize(ldof_group[i]);
725  }
726  // NOTE: in NC meshes, ldof_ltdof generated for the gtopo
727  // groups by ConstructTrueDofs gets overwritten by
728  // BuildParallelConformingInterpolation. Some DOFs that are
729  // seen as true by the conforming code are actually slaves and
730  // end up with a -1 in ldof_ltdof.
731  }
732  }
733 }
734 
736 {
737  GroupCommunicator *gc = new GroupCommunicator(GetGroupTopo());
738  if (NURBSext)
739  {
740  gc->Create(pNURBSext()->ldof_group);
741  }
742  else
743  {
744  GetGroupComm(*gc, 0);
745  }
746  return gc;
747 }
748 
750 {
751  // For non-conforming mesh, synchronization is performed on the cut (aka
752  // "partially conforming") space.
753 
754  MFEM_VERIFY(ldof_marker.Size() == GetVSize(), "invalid in/out array");
755 
756  // implement allreduce(|) as reduce(|) + broadcast
757  gcomm->Reduce<int>(ldof_marker, GroupCommunicator::BitOR);
758  gcomm->Bcast(ldof_marker);
759 }
760 
762  Array<int> &ess_dofs,
763  int component) const
764 {
765  FiniteElementSpace::GetEssentialVDofs(bdr_attr_is_ess, ess_dofs, component);
766 
767  if (Conforming())
768  {
769  // Make sure that processors without boundary elements mark
770  // their boundary dofs (if they have any).
771  Synchronize(ess_dofs);
772  }
773 }
774 
776  &bdr_attr_is_ess,
777  Array<int> &ess_tdof_list,
778  int component)
779 {
780  Array<int> ess_dofs, true_ess_dofs;
781 
782  GetEssentialVDofs(bdr_attr_is_ess, ess_dofs, component);
783  GetRestrictionMatrix()->BooleanMult(ess_dofs, true_ess_dofs);
784 
785 #ifdef MFEM_DEBUG
786  // Verify that in boolean arithmetic: P^T ess_dofs = R ess_dofs.
787  Array<int> true_ess_dofs2(true_ess_dofs.Size());
789  const int *ess_dofs_data = ess_dofs.HostRead();
790  Pt->BooleanMult(1, ess_dofs_data, 0, true_ess_dofs2);
791  delete Pt;
792  int counter = 0;
793  const int *ted = true_ess_dofs.HostRead();
794  for (int i = 0; i < true_ess_dofs.Size(); i++)
795  {
796  if (bool(ted[i]) != bool(true_ess_dofs2[i])) { counter++; }
797  }
798  MFEM_VERIFY(counter == 0, "internal MFEM error: counter = " << counter);
799 #endif
800 
801  MarkerToList(true_ess_dofs, ess_tdof_list);
802 }
803 
805 {
806  if (Nonconforming())
807  {
808  Dof_TrueDof_Matrix(); // make sure P has been built
809 
810  return ldof_ltdof[ldof]; // NOTE: contains -1 for slaves/DOFs we don't own
811  }
812  else
813  {
814  if (GetGroupTopo().IAmMaster(ldof_group[ldof]))
815  {
816  return ldof_ltdof[ldof];
817  }
818  else
819  {
820  return -1;
821  }
822  }
823 }
824 
826 {
827  if (Nonconforming())
828  {
829  MFEM_VERIFY(ldof_ltdof[ldof] >= 0, "ldof " << ldof << " not a true DOF.");
830 
831  return GetMyTDofOffset() + ldof_ltdof[ldof];
832  }
833  else
834  {
835  if (HYPRE_AssumedPartitionCheck())
836  {
837  return ldof_ltdof[ldof] +
838  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(ldof_group[ldof])];
839  }
840  else
841  {
842  return ldof_ltdof[ldof] +
843  tdof_offsets[GetGroupTopo().GetGroupMasterRank(ldof_group[ldof])];
844  }
845  }
846 }
847 
849 {
850  if (Nonconforming())
851  {
852  MFEM_ABORT("Not implemented for NC mesh.");
853  }
854 
855  if (HYPRE_AssumedPartitionCheck())
856  {
858  {
859  return ldof_ltdof[sldof] +
860  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
861  ldof_group[sldof])] / vdim;
862  }
863  else
864  {
865  return (ldof_ltdof[sldof*vdim] +
866  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
867  ldof_group[sldof*vdim])]) / vdim;
868  }
869  }
870 
872  {
873  return ldof_ltdof[sldof] +
874  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
875  ldof_group[sldof])] / vdim;
876  }
877  else
878  {
879  return (ldof_ltdof[sldof*vdim] +
880  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
881  ldof_group[sldof*vdim])]) / vdim;
882  }
883 }
884 
886 {
887  return HYPRE_AssumedPartitionCheck() ? dof_offsets[0] : dof_offsets[MyRank];
888 }
889 
891 {
892  return HYPRE_AssumedPartitionCheck()? tdof_offsets[0] : tdof_offsets[MyRank];
893 }
894 
896 {
897  if (Conforming())
898  {
899  if (Pconf) { return Pconf; }
900 
901  if (NRanks == 1)
902  {
903  Pconf = new IdentityOperator(GetTrueVSize());
904  }
905  else
906  {
908  {
909  Pconf = new ConformingProlongationOperator(*this);
910  }
911  else
912  {
913  Pconf = new DeviceConformingProlongationOperator(*this);
914  }
915  }
916  return Pconf;
917  }
918  else
919  {
920  return Dof_TrueDof_Matrix();
921  }
922 }
923 
925 {
926  if (num_face_nbr_dofs >= 0) { return; }
927 
928  pmesh->ExchangeFaceNbrData();
929 
930  int num_face_nbrs = pmesh->GetNFaceNeighbors();
931 
932  if (num_face_nbrs == 0)
933  {
934  num_face_nbr_dofs = 0;
935  return;
936  }
937 
938  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
939  MPI_Request *send_requests = requests;
940  MPI_Request *recv_requests = requests + num_face_nbrs;
941  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
942 
943  Array<int> ldofs;
944  Array<int> ldof_marker(GetVSize());
945  ldof_marker = -1;
946 
947  Table send_nbr_elem_dof;
948 
949  send_nbr_elem_dof.MakeI(pmesh->send_face_nbr_elements.Size_of_connections());
950  send_face_nbr_ldof.MakeI(num_face_nbrs);
951  face_nbr_ldof.MakeI(num_face_nbrs);
952  int *send_el_off = pmesh->send_face_nbr_elements.GetI();
953  int *recv_el_off = pmesh->face_nbr_elements_offset;
954  for (int fn = 0; fn < num_face_nbrs; fn++)
955  {
956  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
957  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
958 
959  for (int i = 0; i < num_my_elems; i++)
960  {
961  GetElementVDofs(my_elems[i], ldofs);
962  for (int j = 0; j < ldofs.Size(); j++)
963  {
964  int ldof = (ldofs[j] >= 0 ? ldofs[j] : -1-ldofs[j]);
965 
966  if (ldof_marker[ldof] != fn)
967  {
968  ldof_marker[ldof] = fn;
970  }
971  }
972  send_nbr_elem_dof.AddColumnsInRow(send_el_off[fn] + i, ldofs.Size());
973  }
974 
975  int nbr_rank = pmesh->GetFaceNbrRank(fn);
976  int tag = 0;
977  MPI_Isend(&send_face_nbr_ldof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
978  MyComm, &send_requests[fn]);
979 
980  MPI_Irecv(&face_nbr_ldof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
981  MyComm, &recv_requests[fn]);
982  }
983 
984  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
986 
988 
989  MPI_Waitall(num_face_nbrs, send_requests, statuses);
991 
992  // send/receive the I arrays of send_nbr_elem_dof/face_nbr_element_dof,
993  // respectively (they contain the number of dofs for each face-neighbor
994  // element)
995  face_nbr_element_dof.MakeI(recv_el_off[num_face_nbrs]);
996 
997  int *send_I = send_nbr_elem_dof.GetI();
998  int *recv_I = face_nbr_element_dof.GetI();
999  for (int fn = 0; fn < num_face_nbrs; fn++)
1000  {
1001  int nbr_rank = pmesh->GetFaceNbrRank(fn);
1002  int tag = 0;
1003  MPI_Isend(send_I + send_el_off[fn], send_el_off[fn+1] - send_el_off[fn],
1004  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
1005 
1006  MPI_Irecv(recv_I + recv_el_off[fn], recv_el_off[fn+1] - recv_el_off[fn],
1007  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
1008  }
1009 
1010  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1011  send_nbr_elem_dof.MakeJ();
1012 
1013  ldof_marker = -1;
1014 
1015  for (int fn = 0; fn < num_face_nbrs; fn++)
1016  {
1017  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
1018  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
1019 
1020  for (int i = 0; i < num_my_elems; i++)
1021  {
1022  GetElementVDofs(my_elems[i], ldofs);
1023  for (int j = 0; j < ldofs.Size(); j++)
1024  {
1025  int ldof = (ldofs[j] >= 0 ? ldofs[j] : -1-ldofs[j]);
1026 
1027  if (ldof_marker[ldof] != fn)
1028  {
1029  ldof_marker[ldof] = fn;
1030  send_face_nbr_ldof.AddConnection(fn, ldofs[j]);
1031  }
1032  }
1033  send_nbr_elem_dof.AddConnections(
1034  send_el_off[fn] + i, ldofs, ldofs.Size());
1035  }
1036  }
1038  send_nbr_elem_dof.ShiftUpI();
1039 
1040  // convert the ldof indices in send_nbr_elem_dof
1041  int *send_J = send_nbr_elem_dof.GetJ();
1042  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
1043  {
1044  int num_ldofs = send_face_nbr_ldof.RowSize(fn);
1045  int *ldofs = send_face_nbr_ldof.GetRow(fn);
1046  int j_end = send_I[send_el_off[fn+1]];
1047 
1048  for (int i = 0; i < num_ldofs; i++)
1049  {
1050  int ldof = (ldofs[i] >= 0 ? ldofs[i] : -1-ldofs[i]);
1051  ldof_marker[ldof] = i;
1052  }
1053 
1054  for ( ; j < j_end; j++)
1055  {
1056  int ldof = (send_J[j] >= 0 ? send_J[j] : -1-send_J[j]);
1057  send_J[j] = (send_J[j] >= 0 ? ldof_marker[ldof] : -1-ldof_marker[ldof]);
1058  }
1059  }
1060 
1061  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1063 
1064  // send/receive the J arrays of send_nbr_elem_dof/face_nbr_element_dof,
1065  // respectively (they contain the element dofs in enumeration local for
1066  // the face-neighbor pair)
1067  int *recv_J = face_nbr_element_dof.GetJ();
1068  for (int fn = 0; fn < num_face_nbrs; fn++)
1069  {
1070  int nbr_rank = pmesh->GetFaceNbrRank(fn);
1071  int tag = 0;
1072 
1073  MPI_Isend(send_J + send_I[send_el_off[fn]],
1074  send_I[send_el_off[fn+1]] - send_I[send_el_off[fn]],
1075  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
1076 
1077  MPI_Irecv(recv_J + recv_I[recv_el_off[fn]],
1078  recv_I[recv_el_off[fn+1]] - recv_I[recv_el_off[fn]],
1079  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
1080  }
1081 
1082  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1083 
1084  // shift the J array of face_nbr_element_dof
1085  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
1086  {
1087  int shift = face_nbr_ldof.GetI()[fn];
1088  int j_end = recv_I[recv_el_off[fn+1]];
1089 
1090  for ( ; j < j_end; j++)
1091  {
1092  if (recv_J[j] >= 0)
1093  {
1094  recv_J[j] += shift;
1095  }
1096  else
1097  {
1098  recv_J[j] -= shift;
1099  }
1100  }
1101  }
1102 
1103  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1104 
1105  // send/receive the J arrays of send_face_nbr_ldof/face_nbr_ldof,
1106  // respectively
1107  for (int fn = 0; fn < num_face_nbrs; fn++)
1108  {
1109  int nbr_rank = pmesh->GetFaceNbrRank(fn);
1110  int tag = 0;
1111 
1112  MPI_Isend(send_face_nbr_ldof.GetRow(fn),
1114  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
1115 
1116  MPI_Irecv(face_nbr_ldof.GetRow(fn),
1117  face_nbr_ldof.RowSize(fn),
1118  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
1119  }
1120 
1121  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1122  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1123 
1124  // send my_dof_offset (i.e. my_ldof_offset) to face neighbors and receive
1125  // their offset in dof_face_nbr_offsets, used to define face_nbr_glob_dof_map
1127  Array<HYPRE_Int> dof_face_nbr_offsets(num_face_nbrs);
1128  HYPRE_Int my_dof_offset = GetMyDofOffset();
1129  for (int fn = 0; fn < num_face_nbrs; fn++)
1130  {
1131  int nbr_rank = pmesh->GetFaceNbrRank(fn);
1132  int tag = 0;
1133 
1134  MPI_Isend(&my_dof_offset, 1, HYPRE_MPI_INT, nbr_rank, tag,
1135  MyComm, &send_requests[fn]);
1136 
1137  MPI_Irecv(&dof_face_nbr_offsets[fn], 1, HYPRE_MPI_INT, nbr_rank, tag,
1138  MyComm, &recv_requests[fn]);
1139  }
1140 
1141  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
1142 
1143  // set the array face_nbr_glob_dof_map which holds the global ldof indices of
1144  // the face-neighbor dofs
1145  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
1146  {
1147  for (int j_end = face_nbr_ldof.GetI()[fn+1]; j < j_end; j++)
1148  {
1149  int ldof = face_nbr_ldof.GetJ()[j];
1150  if (ldof < 0)
1151  {
1152  ldof = -1-ldof;
1153  }
1154 
1155  face_nbr_glob_dof_map[j] = dof_face_nbr_offsets[fn] + ldof;
1156  }
1157  }
1158 
1159  MPI_Waitall(num_face_nbrs, send_requests, statuses);
1160 
1161  delete [] statuses;
1162  delete [] requests;
1163 }
1164 
1166  int i, Array<int> &vdofs) const
1167 {
1168  face_nbr_element_dof.GetRow(i, vdofs);
1169 }
1170 
1172 {
1173  // Works for NC mesh where 'i' is an index returned by
1174  // ParMesh::GetSharedFace() such that i >= Mesh::GetNumFaces(), i.e. 'i' is
1175  // the index of a ghost face.
1176  MFEM_ASSERT(Nonconforming() && i >= pmesh->GetNumFaces(), "");
1177  int el1, el2, inf1, inf2;
1178  pmesh->GetFaceElements(i, &el1, &el2);
1179  el2 = -1 - el2;
1180  pmesh->GetFaceInfos(i, &inf1, &inf2);
1181  MFEM_ASSERT(0 <= el2 && el2 < face_nbr_element_dof.Size(), "");
1182  const int nd = face_nbr_element_dof.RowSize(el2);
1183  const int *vol_vdofs = face_nbr_element_dof.GetRow(el2);
1184  const Element *face_nbr_el = pmesh->face_nbr_elements[el2];
1185  Geometry::Type geom = face_nbr_el->GetGeometryType();
1186  const int face_dim = Geometry::Dimension[geom]-1;
1187 
1188  fec->SubDofOrder(geom, face_dim, inf2, vdofs);
1189  // Convert local dofs to local vdofs.
1190  Ordering::DofsToVDofs<Ordering::byNODES>(nd/vdim, vdim, vdofs);
1191  // Convert local vdofs to global vdofs.
1192  for (int j = 0; j < vdofs.Size(); j++)
1193  {
1194  const int ldof = vdofs[j];
1195  vdofs[j] = (ldof >= 0) ? vol_vdofs[ldof] : -1-vol_vdofs[-1-ldof];
1196  }
1197 }
1198 
1200 {
1201  const FiniteElement *FE =
1203  pmesh->face_nbr_elements[i]->GetGeometryType());
1204 
1205  if (NURBSext)
1206  {
1207  mfem_error("ParFiniteElementSpace::GetFaceNbrFE"
1208  " does not support NURBS!");
1209  }
1210  return FE;
1211 }
1212 
1214 {
1215  // Works for NC mesh where 'i' is an index returned by
1216  // ParMesh::GetSharedFace() such that i >= Mesh::GetNumFaces(), i.e. 'i' is
1217  // the index of a ghost face.
1218  // Works in tandem with GetFaceNbrFaceVDofs() defined above.
1219 
1220  MFEM_ASSERT(Nonconforming() && !NURBSext, "");
1221  Geometry::Type face_geom = pmesh->GetFaceGeometryType(i);
1222  return fec->FiniteElementForGeometry(face_geom);
1223 }
1224 
1226 {
1227  hypre_ParCSRMatrix *csrP = (hypre_ParCSRMatrix*)(*P);
1228  hypre_ParCSRMatrixOwnsRowStarts(csrP) = 1;
1229  hypre_ParCSRMatrixOwnsColStarts(csrP) = 1;
1230  P -> StealData();
1231  dof_offsets.LoseData();
1232  tdof_offsets.LoseData();
1233 }
1234 
1235 void ParFiniteElementSpace::ConstructTrueDofs()
1236 {
1237  int i, gr, n = GetVSize();
1238  GroupTopology &gt = pmesh->gtopo;
1239  gcomm = new GroupCommunicator(gt);
1240  Table &group_ldof = gcomm->GroupLDofTable();
1241 
1242  GetGroupComm(*gcomm, 1, &ldof_sign);
1243 
1244  // Define ldof_group and mark ldof_ltdof with
1245  // -1 for ldof that is ours
1246  // -2 for ldof that is in a group with another master
1247  ldof_group.SetSize(n);
1248  ldof_ltdof.SetSize(n);
1249  ldof_group = 0;
1250  ldof_ltdof = -1;
1251 
1252  for (gr = 1; gr < group_ldof.Size(); gr++)
1253  {
1254  const int *ldofs = group_ldof.GetRow(gr);
1255  const int nldofs = group_ldof.RowSize(gr);
1256  for (i = 0; i < nldofs; i++)
1257  {
1258  ldof_group[ldofs[i]] = gr;
1259  }
1260 
1261  if (!gt.IAmMaster(gr)) // we are not the master
1262  {
1263  for (i = 0; i < nldofs; i++)
1264  {
1265  ldof_ltdof[ldofs[i]] = -2;
1266  }
1267  }
1268  }
1269 
1270  // count ltdof_size
1271  ltdof_size = 0;
1272  for (i = 0; i < n; i++)
1273  {
1274  if (ldof_ltdof[i] == -1)
1275  {
1276  ldof_ltdof[i] = ltdof_size++;
1277  }
1278  }
1279  gcomm->SetLTDofTable(ldof_ltdof);
1280 
1281  // have the group masters broadcast their ltdofs to the rest of the group
1282  gcomm->Bcast(ldof_ltdof);
1283 }
1284 
1285 void ParFiniteElementSpace::ConstructTrueNURBSDofs()
1286 {
1287  int n = GetVSize();
1288  GroupTopology &gt = pNURBSext()->gtopo;
1289  gcomm = new GroupCommunicator(gt);
1290 
1291  // pNURBSext()->ldof_group is for scalar space!
1292  if (vdim == 1)
1293  {
1294  ldof_group.MakeRef(pNURBSext()->ldof_group);
1295  }
1296  else
1297  {
1298  const int *scalar_ldof_group = pNURBSext()->ldof_group;
1299  ldof_group.SetSize(n);
1300  for (int i = 0; i < n; i++)
1301  {
1302  ldof_group[i] = scalar_ldof_group[VDofToDof(i)];
1303  }
1304  }
1305 
1306  gcomm->Create(ldof_group);
1307 
1308  // ldof_sign.SetSize(n);
1309  // ldof_sign = 1;
1310  ldof_sign.DeleteAll();
1311 
1312  ltdof_size = 0;
1313  ldof_ltdof.SetSize(n);
1314  for (int i = 0; i < n; i++)
1315  {
1316  if (gt.IAmMaster(ldof_group[i]))
1317  {
1318  ldof_ltdof[i] = ltdof_size;
1319  ltdof_size++;
1320  }
1321  else
1322  {
1323  ldof_ltdof[i] = -2;
1324  }
1325  }
1326  gcomm->SetLTDofTable(ldof_ltdof);
1327 
1328  // have the group masters broadcast their ltdofs to the rest of the group
1329  gcomm->Bcast(ldof_ltdof);
1330 }
1331 
1332 void ParFiniteElementSpace::GetGhostVertexDofs(const MeshId &id,
1333  Array<int> &dofs) const
1334 {
1335  int nv = fec->DofForGeometry(Geometry::POINT);
1336  dofs.SetSize(nv);
1337  for (int j = 0; j < nv; j++)
1338  {
1339  dofs[j] = ndofs + nv*id.index + j;
1340  }
1341 }
1342 
1343 void ParFiniteElementSpace::GetGhostEdgeDofs(const MeshId &edge_id,
1344  Array<int> &dofs) const
1345 {
1346  int nv = fec->DofForGeometry(Geometry::POINT);
1348  dofs.SetSize(2*nv + ne);
1349 
1350  int V[2], ghost = pncmesh->GetNVertices();
1351  pmesh->pncmesh->GetEdgeVertices(edge_id, V);
1352 
1353  for (int i = 0; i < 2; i++)
1354  {
1355  int k = (V[i] < ghost) ? V[i]*nv : (ndofs + (V[i] - ghost)*nv);
1356  for (int j = 0; j < nv; j++)
1357  {
1358  dofs[i*nv + j] = k++;
1359  }
1360  }
1361 
1362  int k = ndofs + ngvdofs + (edge_id.index - pncmesh->GetNEdges())*ne;
1363  for (int j = 0; j < ne; j++)
1364  {
1365  dofs[2*nv + j] = k++;
1366  }
1367 }
1368 
1369 void ParFiniteElementSpace::GetGhostFaceDofs(const MeshId &face_id,
1370  Array<int> &dofs) const
1371 {
1372  int nfv, V[4], E[4], Eo[4];
1373  nfv = pmesh->pncmesh->GetFaceVerticesEdges(face_id, V, E, Eo);
1374 
1375  int nv = fec->DofForGeometry(Geometry::POINT);
1377  int nf_tri = fec->DofForGeometry(Geometry::TRIANGLE);
1378  int nf_quad = fec->DofForGeometry(Geometry::SQUARE);
1379  int nf = (nfv == 3) ? nf_tri : nf_quad;
1380 
1381  dofs.SetSize(nfv*(nv + ne) + nf);
1382 
1383  int offset = 0;
1384  for (int i = 0; i < nfv; i++)
1385  {
1386  int ghost = pncmesh->GetNVertices();
1387  int first = (V[i] < ghost) ? V[i]*nv : (ndofs + (V[i] - ghost)*nv);
1388  for (int j = 0; j < nv; j++)
1389  {
1390  dofs[offset++] = first + j;
1391  }
1392  }
1393 
1394  for (int i = 0; i < nfv; i++)
1395  {
1396  int ghost = pncmesh->GetNEdges();
1397  int first = (E[i] < ghost) ? nvdofs + E[i]*ne
1398  /* */ : ndofs + ngvdofs + (E[i] - ghost)*ne;
1399  const int *ind = fec->DofOrderForOrientation(Geometry::SEGMENT, Eo[i]);
1400  for (int j = 0; j < ne; j++)
1401  {
1402  dofs[offset++] = (ind[j] >= 0) ? (first + ind[j])
1403  /* */ : (-1 - (first + (-1 - ind[j])));
1404  }
1405  }
1406 
1407  const int ghost_face_index = face_id.index - pncmesh->GetNFaces();
1408  int first = ndofs + ngvdofs + ngedofs + nf_quad*ghost_face_index;
1409 
1410  for (int j = 0; j < nf; j++)
1411  {
1412  dofs[offset++] = first + j;
1413  }
1414 }
1415 
1416 void ParFiniteElementSpace::GetGhostDofs(int entity, const MeshId &id,
1417  Array<int> &dofs) const
1418 {
1419  // helper to get ghost vertex, ghost edge or ghost face DOFs
1420  switch (entity)
1421  {
1422  case 0: GetGhostVertexDofs(id, dofs); break;
1423  case 1: GetGhostEdgeDofs(id, dofs); break;
1424  case 2: GetGhostFaceDofs(id, dofs); break;
1425  }
1426 }
1427 
1428 void ParFiniteElementSpace::GetBareDofs(int entity, int index,
1429  Array<int> &dofs) const
1430 {
1431  int ned, ghost, first;
1432  switch (entity)
1433  {
1434  case 0:
1436  ghost = pncmesh->GetNVertices();
1437  first = (index < ghost)
1438  ? index*ned // regular vertex
1439  : ndofs + (index - ghost)*ned; // ghost vertex
1440  break;
1441 
1442  case 1:
1444  ghost = pncmesh->GetNEdges();
1445  first = (index < ghost)
1446  ? nvdofs + index*ned // regular edge
1447  : ndofs + ngvdofs + (index - ghost)*ned; // ghost edge
1448  break;
1449 
1450  default:
1451  Geometry::Type geom = pncmesh->GetFaceGeometry(index);
1452  MFEM_ASSERT(geom == Geometry::SQUARE ||
1453  geom == Geometry::TRIANGLE, "");
1454 
1455  ned = fec->DofForGeometry(geom);
1456  ghost = pncmesh->GetNFaces();
1457 
1458  if (index < ghost) // regular face
1459  {
1460  first = nvdofs + nedofs + (fdofs ? fdofs[index] : index*ned);
1461  }
1462  else // ghost face
1463  {
1464  index -= ghost;
1465  int stride = fec->DofForGeometry(Geometry::SQUARE);
1466  first = ndofs + ngvdofs + ngedofs + index*stride;
1467  }
1468  break;
1469  }
1470 
1471  dofs.SetSize(ned);
1472  for (int i = 0; i < ned; i++)
1473  {
1474  dofs[i] = first + i;
1475  }
1476 }
1477 
1478 int ParFiniteElementSpace::PackDof(int entity, int index, int edof) const
1479 {
1480  // DOFs are ordered as follows:
1481  // vertices | edges | faces | internal | ghost vert. | g. edges | g. faces
1482 
1483  int ghost, ned;
1484  switch (entity)
1485  {
1486  case 0:
1487  ghost = pncmesh->GetNVertices();
1489 
1490  return (index < ghost)
1491  ? index*ned + edof // regular vertex
1492  : ndofs + (index - ghost)*ned + edof; // ghost vertex
1493 
1494  case 1:
1495  ghost = pncmesh->GetNEdges();
1497 
1498  return (index < ghost)
1499  ? nvdofs + index*ned + edof // regular edge
1500  : ndofs + ngvdofs + (index - ghost)*ned + edof; // ghost edge
1501 
1502  default:
1503  ghost = pncmesh->GetNFaces();
1504  ned = fec->DofForGeometry(pncmesh->GetFaceGeometry(index));
1505 
1506  if (index < ghost) // regular face
1507  {
1508  return nvdofs + nedofs + (fdofs ? fdofs[index] : index*ned) + edof;
1509  }
1510  else // ghost face
1511  {
1512  index -= ghost;
1513  int stride = fec->DofForGeometry(Geometry::SQUARE);
1514  return ndofs + ngvdofs + ngedofs + index*stride + edof;
1515  }
1516  }
1517 }
1518 
1519 static int bisect(int* array, int size, int value)
1520 {
1521  int* end = array + size;
1522  int* pos = std::upper_bound(array, end, value);
1523  MFEM_VERIFY(pos != end, "value not found");
1524  return pos - array;
1525 }
1526 
1527 /** Dissect a DOF number to obtain the entity type (0=vertex, 1=edge, 2=face),
1528  * entity index and the DOF number within the entity.
1529  */
1530 void ParFiniteElementSpace::UnpackDof(int dof,
1531  int &entity, int &index, int &edof) const
1532 {
1533  MFEM_VERIFY(dof >= 0, "");
1534  if (dof < ndofs)
1535  {
1536  if (dof < nvdofs) // regular vertex
1537  {
1538  int nv = fec->DofForGeometry(Geometry::POINT);
1539  entity = 0, index = dof / nv, edof = dof % nv;
1540  return;
1541  }
1542  dof -= nvdofs;
1543  if (dof < nedofs) // regular edge
1544  {
1546  entity = 1, index = dof / ne, edof = dof % ne;
1547  return;
1548  }
1549  dof -= nedofs;
1550  if (dof < nfdofs) // regular face
1551  {
1552  if (fdofs) // have mixed faces
1553  {
1554  index = bisect(fdofs+1, mesh->GetNFaces(), dof);
1555  edof = dof - fdofs[index];
1556  }
1557  else // uniform faces
1558  {
1559  int nf = fec->DofForGeometry(pncmesh->GetFaceGeometry(0));
1560  index = dof / nf, edof = dof % nf;
1561  }
1562  entity = 2;
1563  return;
1564  }
1565  MFEM_ABORT("Cannot unpack internal DOF");
1566  }
1567  else
1568  {
1569  dof -= ndofs;
1570  if (dof < ngvdofs) // ghost vertex
1571  {
1572  int nv = fec->DofForGeometry(Geometry::POINT);
1573  entity = 0, index = pncmesh->GetNVertices() + dof / nv, edof = dof % nv;
1574  return;
1575  }
1576  dof -= ngvdofs;
1577  if (dof < ngedofs) // ghost edge
1578  {
1580  entity = 1, index = pncmesh->GetNEdges() + dof / ne, edof = dof % ne;
1581  return;
1582  }
1583  dof -= ngedofs;
1584  if (dof < ngfdofs) // ghost face
1585  {
1586  int stride = fec->DofForGeometry(Geometry::SQUARE);
1587  index = pncmesh->GetNFaces() + dof / stride, edof = dof % stride;
1588  entity = 2;
1589  return;
1590  }
1591  MFEM_ABORT("Out of range DOF.");
1592  }
1593 }
1594 
1595 /** Represents an element of the P matrix. The column number is global and
1596  * corresponds to vector dimension 0. The other dimension columns are offset
1597  * by 'stride'.
1598  */
1599 struct PMatrixElement
1600 {
1601  HYPRE_Int column, stride;
1602  double value;
1603 
1604  PMatrixElement(HYPRE_Int col = 0, HYPRE_Int str = 0, double val = 0)
1605  : column(col), stride(str), value(val) {}
1606 
1607  bool operator<(const PMatrixElement &other) const
1608  { return column < other.column; }
1609 
1610  typedef std::vector<PMatrixElement> List;
1611 };
1612 
1613 /** Represents one row of the P matrix, for the construction code below.
1614  * The row is complete: diagonal and offdiagonal elements are not distinguished.
1615  */
1616 struct PMatrixRow
1617 {
1618  PMatrixElement::List elems;
1619 
1620  /// Add other row, times 'coef'.
1621  void AddRow(const PMatrixRow &other, double coef)
1622  {
1623  elems.reserve(elems.size() + other.elems.size());
1624  for (unsigned i = 0; i < other.elems.size(); i++)
1625  {
1626  const PMatrixElement &oei = other.elems[i];
1627  elems.push_back(
1628  PMatrixElement(oei.column, oei.stride, coef * oei.value));
1629  }
1630  }
1631 
1632  /// Remove duplicate columns and sum their values.
1633  void Collapse()
1634  {
1635  if (!elems.size()) { return; }
1636  std::sort(elems.begin(), elems.end());
1637 
1638  int j = 0;
1639  for (unsigned i = 1; i < elems.size(); i++)
1640  {
1641  if (elems[j].column == elems[i].column)
1642  {
1643  elems[j].value += elems[i].value;
1644  }
1645  else
1646  {
1647  elems[++j] = elems[i];
1648  }
1649  }
1650  elems.resize(j+1);
1651  }
1652 
1653  void write(std::ostream &os, double sign) const
1654  {
1655  bin_io::write<int>(os, elems.size());
1656  for (unsigned i = 0; i < elems.size(); i++)
1657  {
1658  const PMatrixElement &e = elems[i];
1659  bin_io::write<HYPRE_Int>(os, e.column);
1660  bin_io::write<int>(os, e.stride); // truncate HYPRE_Int -> int
1661  bin_io::write<double>(os, e.value * sign);
1662  }
1663  }
1664 
1665  void read(std::istream &is, double sign)
1666  {
1667  elems.resize(bin_io::read<int>(is));
1668  for (unsigned i = 0; i < elems.size(); i++)
1669  {
1670  PMatrixElement &e = elems[i];
1671  e.column = bin_io::read<HYPRE_Int>(is);
1672  e.stride = bin_io::read<int>(is);
1673  e.value = bin_io::read<double>(is) * sign;
1674  }
1675  }
1676 };
1677 
1678 /** Represents a message to another processor containing P matrix rows.
1679  * Used by ParFiniteElementSpace::ParallelConformingInterpolation.
1680  */
1681 class NeighborRowMessage : public VarMessage<314>
1682 {
1683 public:
1684  typedef NCMesh::MeshId MeshId;
1685  typedef ParNCMesh::GroupId GroupId;
1686 
1687  struct RowInfo
1688  {
1689  int entity, index, edof;
1690  GroupId group;
1691  PMatrixRow row;
1692 
1693  RowInfo(int ent, int idx, int edof, GroupId grp, const PMatrixRow &row)
1694  : entity(ent), index(idx), edof(edof), group(grp), row(row) {}
1695 
1696  RowInfo(int ent, int idx, int edof, GroupId grp)
1697  : entity(ent), index(idx), edof(edof), group(grp) {}
1698 
1699  typedef std::vector<RowInfo> List;
1700  };
1701 
1702  NeighborRowMessage() : pncmesh(NULL) {}
1703 
1704  void AddRow(int entity, int index, int edof, GroupId group,
1705  const PMatrixRow &row)
1706  {
1707  rows.push_back(RowInfo(entity, index, edof, group, row));
1708  }
1709 
1710  const RowInfo::List& GetRows() const { return rows; }
1711 
1712  void SetNCMesh(ParNCMesh* pnc) { pncmesh = pnc; }
1713  void SetFEC(const FiniteElementCollection* fec) { this->fec = fec; }
1714 
1715  typedef std::map<int, NeighborRowMessage> Map;
1716 
1717 protected:
1718  RowInfo::List rows;
1719 
1720  ParNCMesh *pncmesh;
1721  const FiniteElementCollection* fec;
1722 
1723  virtual void Encode(int rank);
1724  virtual void Decode(int);
1725 };
1726 
1727 
1728 void NeighborRowMessage::Encode(int rank)
1729 {
1730  std::ostringstream stream;
1731 
1732  Array<MeshId> ent_ids[3];
1733  Array<GroupId> group_ids[3];
1734  Array<int> row_idx[3];
1735 
1736  // encode MeshIds and groups
1737  for (unsigned i = 0; i < rows.size(); i++)
1738  {
1739  const RowInfo &ri = rows[i];
1740  const MeshId &id = pncmesh->GetNCList(ri.entity).LookUp(ri.index);
1741  ent_ids[ri.entity].Append(id);
1742  row_idx[ri.entity].Append(i);
1743  group_ids[ri.entity].Append(ri.group);
1744  }
1745 
1746  Array<GroupId> all_group_ids;
1747  all_group_ids.Reserve(rows.size());
1748  for (int i = 0; i < 3; i++)
1749  {
1750  all_group_ids.Append(group_ids[i]);
1751  }
1752 
1753  pncmesh->AdjustMeshIds(ent_ids, rank);
1754  pncmesh->EncodeMeshIds(stream, ent_ids);
1755  pncmesh->EncodeGroups(stream, all_group_ids);
1756 
1757  // write all rows to the stream
1758  for (int ent = 0; ent < 3; ent++)
1759  {
1760  const Array<MeshId> &ids = ent_ids[ent];
1761  for (int i = 0; i < ids.Size(); i++)
1762  {
1763  const MeshId &id = ids[i];
1764  const RowInfo &ri = rows[row_idx[ent][i]];
1765  MFEM_ASSERT(ent == ri.entity, "");
1766 
1767 #ifdef MFEM_DEBUG_PMATRIX
1768  mfem::out << "Rank " << pncmesh->MyRank << " sending to " << rank
1769  << ": ent " << ri.entity << ", index " << ri.index
1770  << ", edof " << ri.edof << " (id " << id.element << "/"
1771  << int(id.local) << ")" << std::endl;
1772 #endif
1773 
1774  // handle orientation and sign change
1775  int edof = ri.edof;
1776  double s = 1.0;
1777  if (ent == 1)
1778  {
1779  int eo = pncmesh->GetEdgeNCOrientation(id);
1780  const int* ind = fec->DofOrderForOrientation(Geometry::SEGMENT, eo);
1781  if ((edof = ind[edof]) < 0)
1782  {
1783  edof = -1 - edof;
1784  s = -1;
1785  }
1786  }
1787 
1788  bin_io::write<int>(stream, edof);
1789  ri.row.write(stream, s);
1790  }
1791  }
1792 
1793  rows.clear();
1794  stream.str().swap(data);
1795 }
1796 
1797 void NeighborRowMessage::Decode(int rank)
1798 {
1799  std::istringstream stream(data);
1800 
1801  Array<MeshId> ent_ids[3];
1802  Array<GroupId> group_ids;
1803 
1804  // decode vertex/edge/face IDs and groups
1805  pncmesh->DecodeMeshIds(stream, ent_ids);
1806  pncmesh->DecodeGroups(stream, group_ids);
1807 
1808  int nrows = ent_ids[0].Size() + ent_ids[1].Size() + ent_ids[2].Size();
1809  MFEM_ASSERT(nrows == group_ids.Size(), "");
1810 
1811  rows.clear();
1812  rows.reserve(nrows);
1813 
1814  // read rows
1815  for (int ent = 0, gi = 0; ent < 3; ent++)
1816  {
1817  const Array<MeshId> &ids = ent_ids[ent];
1818  for (int i = 0; i < ids.Size(); i++)
1819  {
1820  const MeshId &id = ids[i];
1821  int edof = bin_io::read<int>(stream);
1822 
1823  // handle orientation and sign change
1824  const int *ind = NULL;
1825  if (ent == 1)
1826  {
1827  int eo = pncmesh->GetEdgeNCOrientation(id);
1828  ind = fec->DofOrderForOrientation(Geometry::SEGMENT, eo);
1829  }
1830  else if (ent == 2)
1831  {
1832  Geometry::Type geom = pncmesh->GetFaceGeometry(id.index);
1833  int fo = pncmesh->GetFaceOrientation(id.index);
1834  ind = fec->DofOrderForOrientation(geom, fo);
1835  }
1836 
1837  double s = 1.0;
1838  if (ind && (edof = ind[edof]) < 0)
1839  {
1840  edof = -1 - edof;
1841  s = -1.0;
1842  }
1843 
1844  rows.push_back(RowInfo(ent, id.index, edof, group_ids[gi++]));
1845  rows.back().row.read(stream, s);
1846 
1847 #ifdef MFEM_DEBUG_PMATRIX
1848  mfem::out << "Rank " << pncmesh->MyRank << " receiving from " << rank
1849  << ": ent " << rows.back().entity << ", index "
1850  << rows.back().index << ", edof " << rows.back().edof
1851  << std::endl;
1852 #endif
1853  }
1854  }
1855 }
1856 
1857 void
1858 ParFiniteElementSpace::ScheduleSendRow(const PMatrixRow &row, int dof,
1859  GroupId group_id,
1860  NeighborRowMessage::Map &send_msg) const
1861 {
1862  int ent, idx, edof;
1863  UnpackDof(dof, ent, idx, edof);
1864 
1865  const ParNCMesh::CommGroup &group = pncmesh->GetGroup(group_id);
1866  for (unsigned i = 0; i < group.size(); i++)
1867  {
1868  int rank = group[i];
1869  if (rank != MyRank)
1870  {
1871  NeighborRowMessage &msg = send_msg[rank];
1872  msg.AddRow(ent, idx, edof, group_id, row);
1873  msg.SetNCMesh(pncmesh);
1874  msg.SetFEC(fec);
1875 #ifdef MFEM_PMATRIX_STATS
1876  n_rows_sent++;
1877 #endif
1878  }
1879  }
1880 }
1881 
1882 void ParFiniteElementSpace::ForwardRow(const PMatrixRow &row, int dof,
1883  GroupId group_sent_id, GroupId group_id,
1884  NeighborRowMessage::Map &send_msg) const
1885 {
1886  int ent, idx, edof;
1887  UnpackDof(dof, ent, idx, edof);
1888 
1889  const ParNCMesh::CommGroup &group = pncmesh->GetGroup(group_id);
1890  for (unsigned i = 0; i < group.size(); i++)
1891  {
1892  int rank = group[i];
1893  if (rank != MyRank && !pncmesh->GroupContains(group_sent_id, rank))
1894  {
1895  NeighborRowMessage &msg = send_msg[rank];
1896  GroupId invalid = -1; // to prevent forwarding again
1897  msg.AddRow(ent, idx, edof, invalid, row);
1898  msg.SetNCMesh(pncmesh);
1899  msg.SetFEC(fec);
1900 #ifdef MFEM_PMATRIX_STATS
1901  n_rows_fwd++;
1902 #endif
1903 #ifdef MFEM_DEBUG_PMATRIX
1904  mfem::out << "Rank " << pncmesh->GetMyRank() << " forwarding to "
1905  << rank << ": ent " << ent << ", index" << idx
1906  << ", edof " << edof << std::endl;
1907 #endif
1908  }
1909  }
1910 }
1911 
1912 #ifdef MFEM_DEBUG_PMATRIX
1913 void ParFiniteElementSpace
1914 ::DebugDumpDOFs(std::ostream &os,
1915  const SparseMatrix &deps,
1916  const Array<GroupId> &dof_group,
1917  const Array<GroupId> &dof_owner,
1918  const Array<bool> &finalized) const
1919 {
1920  for (int i = 0; i < dof_group.Size(); i++)
1921  {
1922  os << i << ": ";
1923  if (i < (nvdofs + nedofs + nfdofs) || i >= ndofs)
1924  {
1925  int ent, idx, edof;
1926  UnpackDof(i, ent, idx, edof);
1927 
1928  os << edof << " @ ";
1929  if (i > ndofs) { os << "ghost "; }
1930  switch (ent)
1931  {
1932  case 0: os << "vertex "; break;
1933  case 1: os << "edge "; break;
1934  default: os << "face "; break;
1935  }
1936  os << idx << "; ";
1937 
1938  if (i < deps.Height() && deps.RowSize(i))
1939  {
1940  os << "depends on ";
1941  for (int j = 0; j < deps.RowSize(i); j++)
1942  {
1943  os << deps.GetRowColumns(i)[j] << " ("
1944  << deps.GetRowEntries(i)[j] << ")";
1945  if (j < deps.RowSize(i)-1) { os << ", "; }
1946  }
1947  os << "; ";
1948  }
1949  else
1950  {
1951  os << "no deps; ";
1952  }
1953 
1954  os << "group " << dof_group[i] << " (";
1955  const ParNCMesh::CommGroup &g = pncmesh->GetGroup(dof_group[i]);
1956  for (unsigned j = 0; j < g.size(); j++)
1957  {
1958  if (j) { os << ", "; }
1959  os << g[j];
1960  }
1961 
1962  os << "), owner " << dof_owner[i] << " (rank "
1963  << pncmesh->GetGroup(dof_owner[i])[0] << "); "
1964  << (finalized[i] ? "finalized" : "NOT finalized");
1965  }
1966  else
1967  {
1968  os << "internal";
1969  }
1970  os << "\n";
1971  }
1972 }
1973 #endif
1974 
1975 int ParFiniteElementSpace
1976 ::BuildParallelConformingInterpolation(HypreParMatrix **P, SparseMatrix **R,
1977  Array<HYPRE_Int> &dof_offs,
1978  Array<HYPRE_Int> &tdof_offs,
1979  Array<int> *dof_tdof,
1980  bool partial) const
1981 {
1982  bool dg = (nvdofs == 0 && nedofs == 0 && nfdofs == 0);
1983 
1984 #ifdef MFEM_PMATRIX_STATS
1985  n_msgs_sent = n_msgs_recv = 0;
1986  n_rows_sent = n_rows_recv = n_rows_fwd = 0;
1987 #endif
1988 
1989  // *** STEP 1: build master-slave dependency lists ***
1990 
1991  int total_dofs = ndofs + ngdofs;
1992  SparseMatrix deps(ndofs, total_dofs);
1993 
1994  if (!dg && !partial)
1995  {
1996  Array<int> master_dofs, slave_dofs;
1997 
1998  // loop through *all* master edges/faces, constrain their slaves
1999  for (int entity = 0; entity <= 2; entity++)
2000  {
2001  const NCMesh::NCList &list = pncmesh->GetNCList(entity);
2002  if (!list.masters.Size()) { continue; }
2003 
2004  IsoparametricTransformation T;
2005  DenseMatrix I;
2006 
2007  // process masters that we own or that affect our edges/faces
2008  for (int mi = 0; mi < list.masters.Size(); mi++)
2009  {
2010  const NCMesh::Master &mf = list.masters[mi];
2011 
2012  // get master DOFs
2013  pncmesh->IsGhost(entity, mf.index)
2014  ? GetGhostDofs(entity, mf, master_dofs)
2015  : GetEntityDofs(entity, mf.index, master_dofs);
2016 
2017  if (!master_dofs.Size()) { continue; }
2018 
2019  const FiniteElement* fe = fec->FiniteElementForGeometry(mf.Geom());
2020  if (!fe) { continue; }
2021 
2022  switch (mf.Geom())
2023  {
2024  case Geometry::SQUARE: T.SetFE(&QuadrilateralFE); break;
2025  case Geometry::TRIANGLE: T.SetFE(&TriangleFE); break;
2026  case Geometry::SEGMENT: T.SetFE(&SegmentFE); break;
2027  default: MFEM_ABORT("unsupported geometry");
2028  }
2029 
2030  // constrain slaves that exist in our mesh
2031  for (int si = mf.slaves_begin; si < mf.slaves_end; si++)
2032  {
2033  const NCMesh::Slave &sf = list.slaves[si];
2034  if (pncmesh->IsGhost(entity, sf.index)) { continue; }
2035 
2036  GetEntityDofs(entity, sf.index, slave_dofs, mf.Geom());
2037  if (!slave_dofs.Size()) { continue; }
2038 
2039  list.OrientedPointMatrix(sf, T.GetPointMat());
2040  fe->GetLocalInterpolation(T, I);
2041 
2042  // make each slave DOF dependent on all master DOFs
2043  AddDependencies(deps, master_dofs, slave_dofs, I);
2044  }
2045  }
2046  }
2047 
2048  deps.Finalize();
2049  }
2050 
2051  // *** STEP 2: initialize group and owner ID for each DOF ***
2052 
2053  Array<GroupId> dof_group(total_dofs);
2054  Array<GroupId> dof_owner(total_dofs);
2055  dof_group = 0;
2056  dof_owner = 0;
2057 
2058  if (!dg)
2059  {
2060  Array<int> dofs;
2061 
2062  // initialize dof_group[], dof_owner[]
2063  for (int entity = 0; entity <= 2; entity++)
2064  {
2065  const NCMesh::NCList &list = pncmesh->GetNCList(entity);
2066 
2067  int lsize[3] =
2068  { list.conforming.Size(), list.masters.Size(), list.slaves.Size() };
2069 
2070  for (int l = 0; l < 3; l++)
2071  {
2072  for (int i = 0; i < lsize[l]; i++)
2073  {
2074  const MeshId &id =
2075  (l == 0) ? list.conforming[i] :
2076  (l == 1) ? (const MeshId&) list.masters[i]
2077  /* */ : (const MeshId&) list.slaves[i];
2078 
2079  if (id.index < 0) { continue; }
2080 
2081  GroupId owner = pncmesh->GetEntityOwnerId(entity, id.index);
2082  GroupId group = pncmesh->GetEntityGroupId(entity, id.index);
2083 
2084  GetBareDofs(entity, id.index, dofs);
2085 
2086  for (int j = 0; j < dofs.Size(); j++)
2087  {
2088  int dof = dofs[j];
2089  dof_owner[dof] = owner;
2090  dof_group[dof] = group;
2091  }
2092  }
2093  }
2094  }
2095  }
2096 
2097  // *** STEP 3: count true DOFs and calculate P row/column partitions ***
2098 
2099  Array<bool> finalized(total_dofs);
2100  finalized = false;
2101 
2102  // DOFs that stayed independent and are ours are true DOFs
2103  int num_true_dofs = 0;
2104  for (int i = 0; i < ndofs; i++)
2105  {
2106  if (dof_owner[i] == 0 && deps.RowSize(i) == 0)
2107  {
2108  num_true_dofs++;
2109  finalized[i] = true;
2110  }
2111  }
2112 
2113  // calculate global offsets
2114  HYPRE_Int loc_sizes[2] = { ndofs*vdim, num_true_dofs*vdim };
2115  Array<HYPRE_Int>* offsets[2] = { &dof_offs, &tdof_offs };
2116  pmesh->GenerateOffsets(2, loc_sizes, offsets); // calls MPI_Scan, MPI_Bcast
2117 
2118  HYPRE_Int my_tdof_offset =
2119  tdof_offs[HYPRE_AssumedPartitionCheck() ? 0 : MyRank];
2120 
2121  if (R)
2122  {
2123  // initialize the restriction matrix (also parallel but block-diagonal)
2124  *R = new SparseMatrix(num_true_dofs*vdim, ndofs*vdim);
2125  }
2126  if (dof_tdof)
2127  {
2128  dof_tdof->SetSize(ndofs*vdim);
2129  *dof_tdof = -1;
2130  }
2131 
2132  std::vector<PMatrixRow> pmatrix(total_dofs);
2133 
2134  bool bynodes = (ordering == Ordering::byNODES);
2135  int vdim_factor = bynodes ? 1 : vdim;
2136  int dof_stride = bynodes ? ndofs : 1;
2137  int tdof_stride = bynodes ? num_true_dofs : 1;
2138 
2139  // big container for all messages we send (the list is for iterations)
2140  std::list<NeighborRowMessage::Map> send_msg;
2141  send_msg.push_back(NeighborRowMessage::Map());
2142 
2143  // put identity in P and R for true DOFs, set ldof_ltdof
2144  for (int dof = 0, tdof = 0; dof < ndofs; dof++)
2145  {
2146  if (finalized[dof])
2147  {
2148  pmatrix[dof].elems.push_back(
2149  PMatrixElement(my_tdof_offset + vdim_factor*tdof, tdof_stride, 1.));
2150 
2151  // prepare messages to neighbors with identity rows
2152  if (dof_group[dof] != 0)
2153  {
2154  ScheduleSendRow(pmatrix[dof], dof, dof_group[dof], send_msg.back());
2155  }
2156 
2157  for (int vd = 0; vd < vdim; vd++)
2158  {
2159  int vdof = dof*vdim_factor + vd*dof_stride;
2160  int vtdof = tdof*vdim_factor + vd*tdof_stride;
2161 
2162  if (R) { (*R)->Add(vtdof, vdof, 1.0); }
2163  if (dof_tdof) { (*dof_tdof)[vdof] = vtdof; }
2164  }
2165  tdof++;
2166  }
2167  }
2168 
2169  // send identity rows
2170  NeighborRowMessage::IsendAll(send_msg.back(), MyComm);
2171 #ifdef MFEM_PMATRIX_STATS
2172  n_msgs_sent += send_msg.back().size();
2173 #endif
2174 
2175  if (R) { (*R)->Finalize(); }
2176 
2177  // *** STEP 4: main loop ***
2178 
2179  // a single instance (recv_msg) is reused for all incoming messages
2180  NeighborRowMessage recv_msg;
2181  recv_msg.SetNCMesh(pncmesh);
2182  recv_msg.SetFEC(fec);
2183 
2184  int num_finalized = num_true_dofs;
2185  PMatrixRow buffer;
2186  buffer.elems.reserve(1024);
2187 
2188  while (num_finalized < ndofs)
2189  {
2190  // prepare a new round of send buffers
2191  if (send_msg.back().size())
2192  {
2193  send_msg.push_back(NeighborRowMessage::Map());
2194  }
2195 
2196  // check for incoming messages, receive PMatrixRows
2197  int rank, size;
2198  while (NeighborRowMessage::IProbe(rank, size, MyComm))
2199  {
2200  recv_msg.Recv(rank, size, MyComm);
2201 #ifdef MFEM_PMATRIX_STATS
2202  n_msgs_recv++;
2203  n_rows_recv += recv_msg.GetRows().size();
2204 #endif
2205 
2206  const NeighborRowMessage::RowInfo::List &rows = recv_msg.GetRows();
2207  for (unsigned i = 0; i < rows.size(); i++)
2208  {
2209  const NeighborRowMessage::RowInfo &ri = rows[i];
2210  int dof = PackDof(ri.entity, ri.index, ri.edof);
2211  pmatrix[dof] = ri.row;
2212 
2213  if (dof < ndofs && !finalized[dof]) { num_finalized++; }
2214  finalized[dof] = true;
2215 
2216  if (ri.group >= 0 && dof_group[dof] != ri.group)
2217  {
2218  // the sender didn't see the complete group, forward the message
2219  ForwardRow(ri.row, dof, ri.group, dof_group[dof], send_msg.back());
2220  }
2221  }
2222  }
2223 
2224  // finalize all rows that can currently be finalized
2225  bool done = false;
2226  while (!done)
2227  {
2228  done = true;
2229  for (int dof = 0; dof < ndofs; dof++)
2230  {
2231  if (finalized[dof]) { continue; }
2232 
2233  bool owned = (dof_owner[dof] == 0);
2234  bool shared = (dof_group[dof] != 0);
2235 
2236  if (owned && DofFinalizable(dof, finalized, deps))
2237  {
2238  const int* dep_col = deps.GetRowColumns(dof);
2239  const double* dep_coef = deps.GetRowEntries(dof);
2240  int num_dep = deps.RowSize(dof);
2241 
2242  // form linear combination of rows
2243  buffer.elems.clear();
2244  for (int j = 0; j < num_dep; j++)
2245  {
2246  buffer.AddRow(pmatrix[dep_col[j]], dep_coef[j]);
2247  }
2248  buffer.Collapse();
2249  pmatrix[dof] = buffer;
2250 
2251  finalized[dof] = true;
2252  num_finalized++;
2253  done = false;
2254 
2255  // send row to neighbors who need it
2256  if (shared)
2257  {
2258  ScheduleSendRow(pmatrix[dof], dof, dof_group[dof],
2259  send_msg.back());
2260  }
2261  }
2262  }
2263  }
2264 
2265 #ifdef MFEM_DEBUG_PMATRIX
2266  /*static int dump = 0;
2267  if (dump < 10)
2268  {
2269  char fname[100];
2270  sprintf(fname, "dofs%02d.txt", MyRank);
2271  std::ofstream f(fname);
2272  DebugDumpDOFs(f, deps, dof_group, dof_owner, finalized);
2273  dump++;
2274  }*/
2275 #endif
2276 
2277  // send current batch of messages
2278  NeighborRowMessage::IsendAll(send_msg.back(), MyComm);
2279 #ifdef MFEM_PMATRIX_STATS
2280  n_msgs_sent += send_msg.back().size();
2281 #endif
2282  }
2283 
2284  if (P)
2285  {
2286  *P = MakeVDimHypreMatrix(pmatrix, ndofs, num_true_dofs,
2287  dof_offs, tdof_offs);
2288  }
2289 
2290  // clean up possible remaining messages in the queue to avoid receiving
2291  // them erroneously in the next run
2292  int rank, size;
2293  while (NeighborRowMessage::IProbe(rank, size, MyComm))
2294  {
2295  recv_msg.RecvDrop(rank, size, MyComm);
2296  }
2297 
2298  // make sure we can discard all send buffers
2299  for (std::list<NeighborRowMessage::Map>::iterator
2300  it = send_msg.begin(); it != send_msg.end(); ++it)
2301  {
2302  NeighborRowMessage::WaitAllSent(*it);
2303  }
2304 
2305 #ifdef MFEM_PMATRIX_STATS
2306  int n_rounds = send_msg.size();
2307  int glob_rounds, glob_msgs_sent, glob_msgs_recv;
2308  int glob_rows_sent, glob_rows_recv, glob_rows_fwd;
2309 
2310  MPI_Reduce(&n_rounds, &glob_rounds, 1, MPI_INT, MPI_SUM, 0, MyComm);
2311  MPI_Reduce(&n_msgs_sent, &glob_msgs_sent, 1, MPI_INT, MPI_SUM, 0, MyComm);
2312  MPI_Reduce(&n_msgs_recv, &glob_msgs_recv, 1, MPI_INT, MPI_SUM, 0, MyComm);
2313  MPI_Reduce(&n_rows_sent, &glob_rows_sent, 1, MPI_INT, MPI_SUM, 0, MyComm);
2314  MPI_Reduce(&n_rows_recv, &glob_rows_recv, 1, MPI_INT, MPI_SUM, 0, MyComm);
2315  MPI_Reduce(&n_rows_fwd, &glob_rows_fwd, 1, MPI_INT, MPI_SUM, 0, MyComm);
2316 
2317  if (MyRank == 0)
2318  {
2319  mfem::out << "P matrix stats (avg per rank): "
2320  << double(glob_rounds)/NRanks << " rounds, "
2321  << double(glob_msgs_sent)/NRanks << " msgs sent, "
2322  << double(glob_msgs_recv)/NRanks << " msgs recv, "
2323  << double(glob_rows_sent)/NRanks << " rows sent, "
2324  << double(glob_rows_recv)/NRanks << " rows recv, "
2325  << double(glob_rows_fwd)/NRanks << " rows forwarded."
2326  << std::endl;
2327  }
2328 #endif
2329 
2330  return num_true_dofs*vdim;
2331 }
2332 
2333 
2334 HypreParMatrix* ParFiniteElementSpace
2335 ::MakeVDimHypreMatrix(const std::vector<PMatrixRow> &rows,
2336  int local_rows, int local_cols,
2337  Array<HYPRE_Int> &row_starts,
2338  Array<HYPRE_Int> &col_starts) const
2339 {
2340  bool assumed = HYPRE_AssumedPartitionCheck();
2341  bool bynodes = (ordering == Ordering::byNODES);
2342 
2343  HYPRE_Int first_col = col_starts[assumed ? 0 : MyRank];
2344  HYPRE_Int next_col = col_starts[assumed ? 1 : MyRank+1];
2345 
2346  // count nonzeros in diagonal/offdiagonal parts
2347  HYPRE_Int nnz_diag = 0, nnz_offd = 0;
2348  std::map<HYPRE_Int, int> col_map;
2349  for (int i = 0; i < local_rows; i++)
2350  {
2351  for (unsigned j = 0; j < rows[i].elems.size(); j++)
2352  {
2353  const PMatrixElement &elem = rows[i].elems[j];
2354  HYPRE_Int col = elem.column;
2355  if (col >= first_col && col < next_col)
2356  {
2357  nnz_diag += vdim;
2358  }
2359  else
2360  {
2361  nnz_offd += vdim;
2362  for (int vd = 0; vd < vdim; vd++)
2363  {
2364  col_map[col] = -1;
2365  col += elem.stride;
2366  }
2367  }
2368  }
2369  }
2370 
2371  // create offd column mapping
2372  HYPRE_Int *cmap = Memory<HYPRE_Int>(col_map.size());
2373  int offd_col = 0;
2374  for (std::map<HYPRE_Int, int>::iterator
2375  it = col_map.begin(); it != col_map.end(); ++it)
2376  {
2377  cmap[offd_col] = it->first;
2378  it->second = offd_col++;
2379  }
2380 
2381  HYPRE_Int *I_diag = Memory<HYPRE_Int>(vdim*local_rows + 1);
2382  HYPRE_Int *I_offd = Memory<HYPRE_Int>(vdim*local_rows + 1);
2383 
2384  HYPRE_Int *J_diag = Memory<HYPRE_Int>(nnz_diag);
2385  HYPRE_Int *J_offd = Memory<HYPRE_Int>(nnz_offd);
2386 
2387  double *A_diag = Memory<double>(nnz_diag);
2388  double *A_offd = Memory<double>(nnz_offd);
2389 
2390  int vdim1 = bynodes ? vdim : 1;
2391  int vdim2 = bynodes ? 1 : vdim;
2392  int vdim_offset = bynodes ? local_cols : 1;
2393 
2394  // copy the diag/offd elements
2395  nnz_diag = nnz_offd = 0;
2396  int vrow = 0;
2397  for (int vd1 = 0; vd1 < vdim1; vd1++)
2398  {
2399  for (int i = 0; i < local_rows; i++)
2400  {
2401  for (int vd2 = 0; vd2 < vdim2; vd2++)
2402  {
2403  I_diag[vrow] = nnz_diag;
2404  I_offd[vrow++] = nnz_offd;
2405 
2406  int vd = bynodes ? vd1 : vd2;
2407  for (unsigned j = 0; j < rows[i].elems.size(); j++)
2408  {
2409  const PMatrixElement &elem = rows[i].elems[j];
2410  if (elem.column >= first_col && elem.column < next_col)
2411  {
2412  J_diag[nnz_diag] = elem.column + vd*vdim_offset - first_col;
2413  A_diag[nnz_diag++] = elem.value;
2414  }
2415  else
2416  {
2417  J_offd[nnz_offd] = col_map[elem.column + vd*elem.stride];
2418  A_offd[nnz_offd++] = elem.value;
2419  }
2420  }
2421  }
2422  }
2423  }
2424  MFEM_ASSERT(vrow == vdim*local_rows, "");
2425  I_diag[vrow] = nnz_diag;
2426  I_offd[vrow] = nnz_offd;
2427 
2428  return new HypreParMatrix(MyComm,
2429  row_starts.Last(), col_starts.Last(),
2430  row_starts.GetData(), col_starts.GetData(),
2431  I_diag, J_diag, A_diag,
2432  I_offd, J_offd, A_offd,
2433  col_map.size(), cmap);
2434 }
2435 
2436 
2437 static HYPRE_Int* make_i_array(int nrows)
2438 {
2439  HYPRE_Int *I = Memory<HYPRE_Int>(nrows+1);
2440  for (int i = 0; i <= nrows; i++) { I[i] = -1; }
2441  return I;
2442 }
2443 
2444 static HYPRE_Int* make_j_array(HYPRE_Int* I, int nrows)
2445 {
2446  int nnz = 0;
2447  for (int i = 0; i < nrows; i++)
2448  {
2449  if (I[i] >= 0) { nnz++; }
2450  }
2451  HYPRE_Int *J = Memory<HYPRE_Int>(nnz);
2452 
2453  I[nrows] = -1;
2454  for (int i = 0, k = 0; i <= nrows; i++)
2455  {
2456  HYPRE_Int col = I[i];
2457  I[i] = k;
2458  if (col >= 0) { J[k++] = col; }
2459  }
2460  return J;
2461 }
2462 
2463 HypreParMatrix*
2464 ParFiniteElementSpace::RebalanceMatrix(int old_ndofs,
2465  const Table* old_elem_dof)
2466 {
2467  MFEM_VERIFY(Nonconforming(), "Only supported for nonconforming meshes.");
2468  MFEM_VERIFY(old_dof_offsets.Size(), "ParFiniteElementSpace::Update needs to "
2469  "be called before ParFiniteElementSpace::RebalanceMatrix");
2470 
2471  HYPRE_Int old_offset = HYPRE_AssumedPartitionCheck()
2472  ? old_dof_offsets[0] : old_dof_offsets[MyRank];
2473 
2474  // send old DOFs of elements we used to own
2475  ParNCMesh* pncmesh = pmesh->pncmesh;
2476  pncmesh->SendRebalanceDofs(old_ndofs, *old_elem_dof, old_offset, this);
2477 
2478  Array<int> dofs;
2479  int vsize = GetVSize();
2480 
2481  const Array<int> &old_index = pncmesh->GetRebalanceOldIndex();
2482  MFEM_VERIFY(old_index.Size() == pmesh->GetNE(),
2483  "Mesh::Rebalance was not called before "
2484  "ParFiniteElementSpace::RebalanceMatrix");
2485 
2486  // prepare the local (diagonal) part of the matrix
2487  HYPRE_Int* i_diag = make_i_array(vsize);
2488  for (int i = 0; i < pmesh->GetNE(); i++)
2489  {
2490  if (old_index[i] >= 0) // we had this element before
2491  {
2492  const int* old_dofs = old_elem_dof->GetRow(old_index[i]);
2493  GetElementDofs(i, dofs);
2494 
2495  for (int vd = 0; vd < vdim; vd++)
2496  {
2497  for (int j = 0; j < dofs.Size(); j++)
2498  {
2499  int row = DofToVDof(dofs[j], vd);
2500  if (row < 0) { row = -1 - row; }
2501 
2502  int col = DofToVDof(old_dofs[j], vd, old_ndofs);
2503  if (col < 0) { col = -1 - col; }
2504 
2505  i_diag[row] = col;
2506  }
2507  }
2508  }
2509  }
2510  HYPRE_Int* j_diag = make_j_array(i_diag, vsize);
2511 
2512  // receive old DOFs for elements we obtained from others in Rebalance
2513  Array<int> new_elements;
2514  Array<long> old_remote_dofs;
2515  pncmesh->RecvRebalanceDofs(new_elements, old_remote_dofs);
2516 
2517  // create the offdiagonal part of the matrix
2518  HYPRE_Int* i_offd = make_i_array(vsize);
2519  for (int i = 0, pos = 0; i < new_elements.Size(); i++)
2520  {
2521  GetElementDofs(new_elements[i], dofs);
2522  const long* old_dofs = &old_remote_dofs[pos];
2523  pos += dofs.Size() * vdim;
2524 
2525  for (int vd = 0; vd < vdim; vd++)
2526  {
2527  for (int j = 0; j < dofs.Size(); j++)
2528  {
2529  int row = DofToVDof(dofs[j], vd);
2530  if (row < 0) { row = -1 - row; }
2531 
2532  if (i_diag[row] == i_diag[row+1]) // diag row empty?
2533  {
2534  i_offd[row] = old_dofs[j + vd * dofs.Size()];
2535  }
2536  }
2537  }
2538  }
2539  HYPRE_Int* j_offd = make_j_array(i_offd, vsize);
2540 
2541  // create the offd column map
2542  int offd_cols = i_offd[vsize];
2543  Array<Pair<HYPRE_Int, int> > cmap_offd(offd_cols);
2544  for (int i = 0; i < offd_cols; i++)
2545  {
2546  cmap_offd[i].one = j_offd[i];
2547  cmap_offd[i].two = i;
2548  }
2549  SortPairs<HYPRE_Int, int>(cmap_offd, offd_cols);
2550 
2551  HYPRE_Int* cmap = Memory<HYPRE_Int>(offd_cols);
2552  for (int i = 0; i < offd_cols; i++)
2553  {
2554  cmap[i] = cmap_offd[i].one;
2555  j_offd[cmap_offd[i].two] = i;
2556  }
2557 
2558  HypreParMatrix *M;
2559  M = new HypreParMatrix(MyComm, MyRank, NRanks, dof_offsets, old_dof_offsets,
2560  i_diag, j_diag, i_offd, j_offd, cmap, offd_cols);
2561  return M;
2562 }
2563 
2564 
2565 struct DerefDofMessage
2566 {
2567  std::vector<HYPRE_Int> dofs;
2568  MPI_Request request;
2569 };
2570 
2571 HypreParMatrix*
2572 ParFiniteElementSpace::ParallelDerefinementMatrix(int old_ndofs,
2573  const Table* old_elem_dof)
2574 {
2575  int nrk = HYPRE_AssumedPartitionCheck() ? 2 : NRanks;
2576 
2577  MFEM_VERIFY(Nonconforming(), "Not implemented for conforming meshes.");
2578  MFEM_VERIFY(old_dof_offsets[nrk], "Missing previous (finer) space.");
2579 
2580 #if 0 // check no longer seems to work with NC tet refinement
2581  MFEM_VERIFY(dof_offsets[nrk] <= old_dof_offsets[nrk],
2582  "Previous space is not finer.");
2583 #endif
2584 
2585  // Note to the reader: please make sure you first read
2586  // FiniteElementSpace::RefinementMatrix, then
2587  // FiniteElementSpace::DerefinementMatrix, and only then this function.
2588  // You have been warned! :-)
2589 
2590  Mesh::GeometryList elem_geoms(*mesh);
2591 
2592  Array<int> dofs, old_dofs, old_vdofs;
2593  Vector row;
2594 
2595  ParNCMesh* pncmesh = pmesh->pncmesh;
2596 
2597  int ldof[Geometry::NumGeom];
2598  for (int i = 0; i < Geometry::NumGeom; i++)
2599  {
2600  ldof[i] = 0;
2601  }
2602  for (int i = 0; i < elem_geoms.Size(); i++)
2603  {
2604  Geometry::Type geom = elem_geoms[i];
2605  ldof[geom] = fec->FiniteElementForGeometry(geom)->GetDof();
2606  }
2607 
2608  const CoarseFineTransformations &dtrans = pncmesh->GetDerefinementTransforms();
2609  const Array<int> &old_ranks = pncmesh->GetDerefineOldRanks();
2610 
2611  std::map<int, DerefDofMessage> messages;
2612 
2613  HYPRE_Int old_offset = HYPRE_AssumedPartitionCheck()
2614  ? old_dof_offsets[0] : old_dof_offsets[MyRank];
2615 
2616  // communicate DOFs for derefinements that straddle processor boundaries,
2617  // note that this is infrequent due to the way elements are ordered
2618  for (int k = 0; k < dtrans.embeddings.Size(); k++)
2619  {
2620  const Embedding &emb = dtrans.embeddings[k];
2621 
2622  int fine_rank = old_ranks[k];
2623  int coarse_rank = (emb.parent < 0) ? (-1 - emb.parent)
2624  : pncmesh->ElementRank(emb.parent);
2625 
2626  if (coarse_rank != MyRank && fine_rank == MyRank)
2627  {
2628  old_elem_dof->GetRow(k, dofs);
2629  DofsToVDofs(dofs, old_ndofs);
2630 
2631  DerefDofMessage &msg = messages[k];
2632  msg.dofs.resize(dofs.Size());
2633  for (int i = 0; i < dofs.Size(); i++)
2634  {
2635  msg.dofs[i] = old_offset + dofs[i];
2636  }
2637 
2638  MPI_Isend(&msg.dofs[0], msg.dofs.size(), HYPRE_MPI_INT,
2639  coarse_rank, 291, MyComm, &msg.request);
2640  }
2641  else if (coarse_rank == MyRank && fine_rank != MyRank)
2642  {
2643  MFEM_ASSERT(emb.parent >= 0, "");
2644  Geometry::Type geom = mesh->GetElementBaseGeometry(emb.parent);
2645 
2646  DerefDofMessage &msg = messages[k];
2647  msg.dofs.resize(ldof[geom]*vdim);
2648 
2649  MPI_Irecv(&msg.dofs[0], ldof[geom]*vdim, HYPRE_MPI_INT,
2650  fine_rank, 291, MyComm, &msg.request);
2651  }
2652  // TODO: coalesce Isends/Irecvs to the same rank. Typically, on uniform
2653  // derefinement, there should be just one send to MyRank-1 and one recv
2654  // from MyRank+1
2655  }
2656 
2657  DenseTensor localR[Geometry::NumGeom];
2658  for (int i = 0; i < elem_geoms.Size(); i++)
2659  {
2660  GetLocalDerefinementMatrices(elem_geoms[i], localR[elem_geoms[i]]);
2661  }
2662 
2663  // create the diagonal part of the derefinement matrix
2664  SparseMatrix *diag = new SparseMatrix(ndofs*vdim, old_ndofs*vdim);
2665 
2666  Array<char> mark(diag->Height());
2667  mark = 0;
2668 
2669  for (int k = 0; k < dtrans.embeddings.Size(); k++)
2670  {
2671  const Embedding &emb = dtrans.embeddings[k];
2672  if (emb.parent < 0) { continue; }
2673 
2674  int coarse_rank = pncmesh->ElementRank(emb.parent);
2675  int fine_rank = old_ranks[k];
2676 
2677  if (coarse_rank == MyRank && fine_rank == MyRank)
2678  {
2679  Geometry::Type geom = mesh->GetElementBaseGeometry(emb.parent);
2680  DenseMatrix &lR = localR[geom](emb.matrix);
2681 
2682  elem_dof->GetRow(emb.parent, dofs);
2683  old_elem_dof->GetRow(k, old_dofs);
2684 
2685  for (int vd = 0; vd < vdim; vd++)
2686  {
2687  old_dofs.Copy(old_vdofs);
2688  DofsToVDofs(vd, old_vdofs, old_ndofs);
2689 
2690  for (int i = 0; i < lR.Height(); i++)
2691  {
2692  if (!std::isfinite(lR(i, 0))) { continue; }
2693 
2694  int r = DofToVDof(dofs[i], vd);
2695  int m = (r >= 0) ? r : (-1 - r);
2696 
2697  if (!mark[m])
2698  {
2699  lR.GetRow(i, row);
2700  diag->SetRow(r, old_vdofs, row);
2701  mark[m] = 1;
2702  }
2703  }
2704  }
2705  }
2706  }
2707  diag->Finalize();
2708 
2709  // wait for all sends/receives to complete
2710  for (auto it = messages.begin(); it != messages.end(); ++it)
2711  {
2712  MPI_Wait(&it->second.request, MPI_STATUS_IGNORE);
2713  }
2714 
2715  // create the offdiagonal part of the derefinement matrix
2716  SparseMatrix *offd = new SparseMatrix(ndofs*vdim, 1);
2717 
2718  std::map<HYPRE_Int, int> col_map;
2719  for (int k = 0; k < dtrans.embeddings.Size(); k++)
2720  {
2721  const Embedding &emb = dtrans.embeddings[k];
2722  if (emb.parent < 0) { continue; }
2723 
2724  int coarse_rank = pncmesh->ElementRank(emb.parent);
2725  int fine_rank = old_ranks[k];
2726 
2727  if (coarse_rank == MyRank && fine_rank != MyRank)
2728  {
2729  Geometry::Type geom = mesh->GetElementBaseGeometry(emb.parent);
2730  DenseMatrix &lR = localR[geom](emb.matrix);
2731 
2732  elem_dof->GetRow(emb.parent, dofs);
2733 
2734  DerefDofMessage &msg = messages[k];
2735  MFEM_ASSERT(msg.dofs.size(), "");
2736 
2737  for (int vd = 0; vd < vdim; vd++)
2738  {
2739  MFEM_ASSERT(ldof[geom], "");
2740  HYPRE_Int* remote_dofs = &msg.dofs[vd*ldof[geom]];
2741 
2742  for (int i = 0; i < lR.Height(); i++)
2743  {
2744  if (!std::isfinite(lR(i, 0))) { continue; }
2745 
2746  int r = DofToVDof(dofs[i], vd);
2747  int m = (r >= 0) ? r : (-1 - r);
2748 
2749  if (!mark[m])
2750  {
2751  lR.GetRow(i, row);
2752  MFEM_ASSERT(ldof[geom] == row.Size(), "");
2753  for (int j = 0; j < ldof[geom]; j++)
2754  {
2755  if (row[j] == 0.0) { continue; } // NOTE: lR thresholded
2756  int &lcol = col_map[remote_dofs[j]];
2757  if (!lcol) { lcol = col_map.size(); }
2758  offd->_Set_(m, lcol-1, row[j]);
2759  }
2760  mark[m] = 1;
2761  }
2762  }
2763  }
2764  }
2765  }
2766  messages.clear();
2767  offd->Finalize(0);
2768  offd->SetWidth(col_map.size());
2769 
2770  // create offd column mapping for use by hypre
2771  HYPRE_Int *cmap = Memory<HYPRE_Int>(offd->Width());
2772  for (std::map<HYPRE_Int, int>::iterator
2773  it = col_map.begin(); it != col_map.end(); ++it)
2774  {
2775  cmap[it->second-1] = it->first;
2776  }
2777 
2778  // reorder offd columns so that 'cmap' is monotonic
2779  // NOTE: this is easier and probably faster (offd is small) than making
2780  // sure cmap is determined and sorted before the offd matrix is created
2781  {
2782  int width = offd->Width();
2783  Array<Pair<HYPRE_Int, int> > reorder(width);
2784  for (int i = 0; i < width; i++)
2785  {
2786  reorder[i].one = cmap[i];
2787  reorder[i].two = i;
2788  }
2789  reorder.Sort();
2790 
2791  Array<int> reindex(width);
2792  for (int i = 0; i < width; i++)
2793  {
2794  reindex[reorder[i].two] = i;
2795  cmap[i] = reorder[i].one;
2796  }
2797 
2798  int *J = offd->GetJ();
2799  for (int i = 0; i < offd->NumNonZeroElems(); i++)
2800  {
2801  J[i] = reindex[J[i]];
2802  }
2803  offd->SortColumnIndices();
2804  }
2805 
2806  HypreParMatrix* R;
2807  R = new HypreParMatrix(MyComm, dof_offsets[nrk], old_dof_offsets[nrk],
2808  dof_offsets, old_dof_offsets, diag, offd, cmap);
2809 
2810 #ifndef HYPRE_BIGINT
2811  diag->LoseData();
2812  offd->LoseData();
2813 #else
2814  diag->SetDataOwner(false);
2815  offd->SetDataOwner(false);
2816 #endif
2817  delete diag;
2818  delete offd;
2819 
2820  R->SetOwnerFlags(3, 3, 1);
2821 
2822  return R;
2823 }
2824 
2825 void ParFiniteElementSpace::Destroy()
2826 {
2827  ldof_group.DeleteAll();
2828  ldof_ltdof.DeleteAll();
2829  dof_offsets.DeleteAll();
2830  tdof_offsets.DeleteAll();
2831  tdof_nb_offsets.DeleteAll();
2832  // preserve old_dof_offsets
2833  ldof_sign.DeleteAll();
2834 
2835  delete P; P = NULL;
2836  delete Pconf; Pconf = NULL;
2837  delete R; R = NULL;
2838 
2839  delete gcomm; gcomm = NULL;
2840 
2841  num_face_nbr_dofs = -1;
2843  face_nbr_ldof.Clear();
2846 }
2847 
2849  const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
2850 {
2853  GetTransferOperator(coarse_fes, Tgf);
2854  Dof_TrueDof_Matrix(); // Make sure R is built - we need R in all cases.
2855  if (T.Type() == Operator::Hypre_ParCSR)
2856  {
2857  const ParFiniteElementSpace *c_pfes =
2858  dynamic_cast<const ParFiniteElementSpace *>(&coarse_fes);
2859  MFEM_ASSERT(c_pfes != NULL, "coarse_fes must be a parallel space");
2860  SparseMatrix *RA = mfem::Mult(*R, *Tgf.As<SparseMatrix>());
2861  Tgf.Clear();
2862  T.Reset(c_pfes->Dof_TrueDof_Matrix()->
2863  LeftDiagMult(*RA, GetTrueDofOffsets()));
2864  delete RA;
2865  }
2866  else
2867  {
2868  T.Reset(new TripleProductOperator(R, Tgf.Ptr(),
2869  coarse_fes.GetProlongationMatrix(),
2870  false, Tgf.OwnsOperator(), false));
2871  Tgf.SetOperatorOwner(false);
2872  }
2873 }
2874 
2875 void ParFiniteElementSpace::Update(bool want_transform)
2876 {
2877  if (mesh->GetSequence() == sequence)
2878  {
2879  return; // no need to update, no-op
2880  }
2881  if (want_transform && mesh->GetSequence() != sequence + 1)
2882  {
2883  MFEM_ABORT("Error in update sequence. Space needs to be updated after "
2884  "each mesh modification.");
2885  }
2886  sequence = mesh->GetSequence();
2887 
2888  if (NURBSext)
2889  {
2890  UpdateNURBS();
2891  return;
2892  }
2893 
2894  Table* old_elem_dof = NULL;
2895  int old_ndofs;
2896 
2897  // save old DOF table
2898  if (want_transform)
2899  {
2900  old_elem_dof = elem_dof;
2901  elem_dof = NULL;
2902  old_ndofs = ndofs;
2903  Swap(dof_offsets, old_dof_offsets);
2904  }
2905 
2906  Destroy();
2907  FiniteElementSpace::Destroy(); // calls Th.Clear()
2908 
2910  Construct();
2911 
2913 
2914  if (want_transform)
2915  {
2916  // calculate appropriate GridFunction transformation
2917  switch (mesh->GetLastOperation())
2918  {
2919  case Mesh::REFINE:
2920  {
2921  if (Th.Type() != Operator::MFEM_SPARSEMAT)
2922  {
2923  Th.Reset(new RefinementOperator(this, old_elem_dof, old_ndofs));
2924  // The RefinementOperator takes ownership of 'old_elem_dofs', so
2925  // we no longer own it:
2926  old_elem_dof = NULL;
2927  }
2928  else
2929  {
2930  // calculate fully assembled matrix
2931  Th.Reset(RefinementMatrix(old_ndofs, old_elem_dof));
2932  }
2933  break;
2934  }
2935 
2936  case Mesh::DEREFINE:
2937  {
2938  Th.Reset(ParallelDerefinementMatrix(old_ndofs, old_elem_dof));
2939  if (Nonconforming())
2940  {
2941  Th.SetOperatorOwner(false);
2942  Th.Reset(new TripleProductOperator(P, R, Th.Ptr(),
2943  false, false, true));
2944  }
2945  break;
2946  }
2947 
2948  case Mesh::REBALANCE:
2949  {
2950  Th.Reset(RebalanceMatrix(old_ndofs, old_elem_dof));
2951  break;
2952  }
2953 
2954  default:
2955  break;
2956  }
2957 
2958  delete old_elem_dof;
2959  }
2960 }
2961 
2962 
2964  const ParFiniteElementSpace &pfes)
2965  : Operator(pfes.GetVSize(), pfes.GetTrueVSize()),
2966  external_ldofs(),
2967  gc(pfes.GroupComm())
2968 {
2969  MFEM_VERIFY(pfes.Conforming(), "");
2970  const Table &group_ldof = gc.GroupLDofTable();
2972  for (int gr = 1; gr < group_ldof.Size(); gr++)
2973  {
2974  if (!gc.GetGroupTopology().IAmMaster(gr))
2975  {
2976  external_ldofs.Append(group_ldof.GetRow(gr), group_ldof.RowSize(gr));
2977  }
2978  }
2979  external_ldofs.Sort();
2980  MFEM_ASSERT(external_ldofs.Size() == Height()-Width(), "");
2981 #ifdef MFEM_DEBUG
2982  for (int j = 1; j < external_ldofs.Size(); j++)
2983  {
2984  // Check for repeated ldofs.
2985  MFEM_VERIFY(external_ldofs[j-1] < external_ldofs[j], "");
2986  }
2987  int j = 0;
2988  for (int i = 0; i < external_ldofs.Size(); i++)
2989  {
2990  const int end = external_ldofs[i];
2991  for ( ; j < end; j++)
2992  {
2993  MFEM_VERIFY(j-i == pfes.GetLocalTDofNumber(j), "");
2994  }
2995  j = end+1;
2996  }
2997  for ( ; j < Height(); j++)
2998  {
2999  MFEM_VERIFY(j-external_ldofs.Size() == pfes.GetLocalTDofNumber(j), "");
3000  }
3001  // gc.PrintInfo();
3002  // pfes.Dof_TrueDof_Matrix()->PrintCommPkg();
3003 #endif
3004 }
3005 
3007 {
3008  MFEM_ASSERT(x.Size() == Width(), "");
3009  MFEM_ASSERT(y.Size() == Height(), "");
3010 
3011  const double *xdata = x.HostRead();
3012  double *ydata = y.HostWrite();
3013  const int m = external_ldofs.Size();
3014 
3015  const int in_layout = 2; // 2 - input is ltdofs array
3016  gc.BcastBegin(const_cast<double*>(xdata), in_layout);
3017 
3018  int j = 0;
3019  for (int i = 0; i < m; i++)
3020  {
3021  const int end = external_ldofs[i];
3022  std::copy(xdata+j-i, xdata+end-i, ydata+j);
3023  j = end+1;
3024  }
3025  std::copy(xdata+j-m, xdata+Width(), ydata+j);
3026 
3027  const int out_layout = 0; // 0 - output is ldofs array
3028  gc.BcastEnd(ydata, out_layout);
3029 }
3030 
3032  const Vector &x, Vector &y) const
3033 {
3034  MFEM_ASSERT(x.Size() == Height(), "");
3035  MFEM_ASSERT(y.Size() == Width(), "");
3036 
3037  const double *xdata = x.HostRead();
3038  double *ydata = y.HostWrite();
3039  const int m = external_ldofs.Size();
3040 
3041  gc.ReduceBegin(xdata);
3042 
3043  int j = 0;
3044  for (int i = 0; i < m; i++)
3045  {
3046  const int end = external_ldofs[i];
3047  std::copy(xdata+j, xdata+end, ydata+j-i);
3048  j = end+1;
3049  }
3050  std::copy(xdata+j, xdata+Height(), ydata+j-m);
3051 
3052  const int out_layout = 2; // 2 - output is an array on all ltdofs
3053  gc.ReduceEnd<double>(ydata, out_layout, GroupCommunicator::Sum);
3054 }
3055 
3057  const ParFiniteElementSpace &pfes) :
3059  mpi_gpu_aware(Device::GetGPUAwareMPI())
3060 {
3061  MFEM_ASSERT(pfes.Conforming(), "internal error");
3062  const SparseMatrix *R = pfes.GetRestrictionMatrix();
3063  MFEM_ASSERT(R->Finalized(), "");
3064  const int tdofs = R->Height();
3065  MFEM_ASSERT(tdofs == pfes.GetTrueVSize(), "");
3066  MFEM_ASSERT(tdofs == R->HostReadI()[tdofs], "");
3067  ltdof_ldof = Array<int>(const_cast<int*>(R->HostReadJ()), tdofs);
3068  ltdof_ldof.UseDevice();
3069  {
3070  Table nbr_ltdof;
3071  gc.GetNeighborLTDofTable(nbr_ltdof);
3072  const int nb_connections = nbr_ltdof.Size_of_connections();
3073  shr_ltdof.SetSize(nb_connections);
3074  shr_ltdof.CopyFrom(nbr_ltdof.GetJ());
3075  shr_buf.SetSize(nb_connections);
3076  shr_buf.UseDevice(true);
3077  shr_buf_offsets = nbr_ltdof.GetIMemory();
3078  {
3079  Array<int> shr_ltdof(nbr_ltdof.GetJ(), nb_connections);
3080  Array<int> unique_ltdof(shr_ltdof);
3081  unique_ltdof.Sort();
3082  unique_ltdof.Unique();
3083  // Note: the next loop modifies the J array of nbr_ltdof
3084  for (int i = 0; i < shr_ltdof.Size(); i++)
3085  {
3086  shr_ltdof[i] = unique_ltdof.FindSorted(shr_ltdof[i]);
3087  MFEM_ASSERT(shr_ltdof[i] != -1, "internal error");
3088  }
3089  Table unique_shr;
3090  Transpose(shr_ltdof, unique_shr, unique_ltdof.Size());
3091  unq_ltdof = Array<int>(unique_ltdof, unique_ltdof.Size());
3092  unq_shr_i = Array<int>(unique_shr.GetI(), unique_shr.Size()+1);
3093  unq_shr_j = Array<int>(unique_shr.GetJ(), unique_shr.Size_of_connections());
3094  }
3095  nbr_ltdof.GetJMemory().Delete();
3096  nbr_ltdof.LoseData();
3097  }
3098  {
3099  Table nbr_ldof;
3100  gc.GetNeighborLDofTable(nbr_ldof);
3101  const int nb_connections = nbr_ldof.Size_of_connections();
3102  ext_ldof.SetSize(nb_connections);
3103  ext_ldof.CopyFrom(nbr_ldof.GetJ());
3104  ext_buf.SetSize(nb_connections);
3105  ext_buf.UseDevice(true);
3106  ext_buf_offsets = nbr_ldof.GetIMemory();
3107  nbr_ldof.GetJMemory().Delete();
3108  nbr_ldof.LoseData();
3109  }
3110  const GroupTopology &gtopo = gc.GetGroupTopology();
3111  int req_counter = 0;
3112  for (int nbr = 1; nbr < gtopo.GetNumNeighbors(); nbr++)
3113  {
3114  const int send_offset = shr_buf_offsets[nbr];
3115  const int send_size = shr_buf_offsets[nbr+1] - send_offset;
3116  if (send_size > 0) { req_counter++; }
3117 
3118  const int recv_offset = ext_buf_offsets[nbr];
3119  const int recv_size = ext_buf_offsets[nbr+1] - recv_offset;
3120  if (recv_size > 0) { req_counter++; }
3121  }
3122  requests = new MPI_Request[req_counter];
3123 }
3124 
3125 static void ExtractSubVector(const int N,
3126  const Array<int> &indices,
3127  const Vector &in, Vector &out)
3128 {
3129  auto y = out.Write();
3130  const auto x = in.Read();
3131  const auto I = indices.Read();
3132  MFEM_FORALL(i, N, y[i] = x[I[i]];); // indices can be repeated
3133 }
3134 
3136  const Vector &x) const
3137 {
3138  // shr_buf[i] = src[shr_ltdof[i]]
3139  if (shr_ltdof.Size() == 0) { return; }
3140  ExtractSubVector(shr_ltdof.Size(), shr_ltdof, x, shr_buf);
3141  // If the above kernel is executed asynchronously, we should wait for it to
3142  // complete
3143  if (mpi_gpu_aware) { MFEM_STREAM_SYNC; }
3144 }
3145 
3146 static void SetSubVector(const int N,
3147  const Array<int> &indices,
3148  const Vector &in, Vector &out)
3149 {
3150  auto y = out.Write();
3151  const auto x = in.Read();
3152  const auto I = indices.Read();
3153  MFEM_FORALL(i, N, y[I[i]] = x[i];);
3154 }
3155 
3157  const Vector &x, Vector &y) const
3158 {
3159  // dst[ltdof_ldof[i]] = src[i]
3160  if (ltdof_ldof.Size() == 0) { return; }
3161  SetSubVector(ltdof_ldof.Size(), ltdof_ldof, x, y);
3162 }
3163 
3165  Vector &y) const
3166 {
3167  // dst[ext_ldof[i]] = ext_buf[i]
3168  if (ext_ldof.Size() == 0) { return; }
3169  SetSubVector(ext_ldof.Size(), ext_ldof, ext_buf, y);
3170 }
3171 
3173  Vector &y) const
3174 {
3175  const GroupTopology &gtopo = gc.GetGroupTopology();
3176  BcastBeginCopy(x); // copy to 'shr_buf'
3177  int req_counter = 0;
3178  for (int nbr = 1; nbr < gtopo.GetNumNeighbors(); nbr++)
3179  {
3180  const int send_offset = shr_buf_offsets[nbr];
3181  const int send_size = shr_buf_offsets[nbr+1] - send_offset;
3182  if (send_size > 0)
3183  {
3184  auto send_buf = mpi_gpu_aware ? shr_buf.Read() : shr_buf.HostRead();
3185  MPI_Isend(send_buf + send_offset, send_size, MPI_DOUBLE,
3186  gtopo.GetNeighborRank(nbr), 41822,
3187  gtopo.GetComm(), &requests[req_counter++]);
3188  }
3189  const int recv_offset = ext_buf_offsets[nbr];
3190  const int recv_size = ext_buf_offsets[nbr+1] - recv_offset;
3191  if (recv_size > 0)
3192  {
3193  auto recv_buf = mpi_gpu_aware ? ext_buf.Write() : ext_buf.HostWrite();
3194  MPI_Irecv(recv_buf + recv_offset, recv_size, MPI_DOUBLE,
3195  gtopo.GetNeighborRank(nbr), 41822,
3196  gtopo.GetComm(), &requests[req_counter++]);
3197  }
3198  }
3199  BcastLocalCopy(x, y);
3200  MPI_Waitall(req_counter, requests, MPI_STATUSES_IGNORE);
3201  BcastEndCopy(y); // copy from 'ext_buf'
3202 }
3203 
3205 {
3206  delete [] requests;
3209 }
3210 
3212  const Vector &x) const
3213 {
3214  // ext_buf[i] = src[ext_ldof[i]]
3215  if (ext_ldof.Size() == 0) { return; }
3216  ExtractSubVector(ext_ldof.Size(), ext_ldof, x, ext_buf);
3217  // If the above kernel is executed asynchronously, we should wait for it to
3218  // complete
3219  if (mpi_gpu_aware) { MFEM_STREAM_SYNC; }
3220 }
3221 
3223  const Vector &x, Vector &y) const
3224 {
3225  // dst[i] = src[ltdof_ldof[i]]
3226  if (ltdof_ldof.Size() == 0) { return; }
3227  ExtractSubVector(ltdof_ldof.Size(), ltdof_ldof, x, y);
3228 }
3229 
3230 static void AddSubVector(const int num_unique_dst_indices,
3231  const Array<int> &unique_dst_indices,
3232  const Array<int> &unique_to_src_offsets,
3233  const Array<int> &unique_to_src_indices,
3234  const Vector &src,
3235  Vector &dst)
3236 {
3237  auto y = dst.Write();
3238  const auto x = src.Read();
3239  const auto DST_I = unique_dst_indices.Read();
3240  const auto SRC_O = unique_to_src_offsets.Read();
3241  const auto SRC_I = unique_to_src_indices.Read();
3242  MFEM_FORALL(i, num_unique_dst_indices,
3243  {
3244  const int dst_idx = DST_I[i];
3245  double sum = y[dst_idx];
3246  const int end = SRC_O[i+1];
3247  for (int j = SRC_O[i]; j != end; ++j) { sum += x[SRC_I[j]]; }
3248  y[dst_idx] = sum;
3249  });
3250 }
3251 
3253 {
3254  // dst[shr_ltdof[i]] += shr_buf[i]
3255  const int unq_ltdof_size = unq_ltdof.Size();
3256  if (unq_ltdof_size == 0) { return; }
3257  AddSubVector(unq_ltdof_size, unq_ltdof, unq_shr_i, unq_shr_j, shr_buf, y);
3258 }
3259 
3261  Vector &y) const
3262 {
3263  const GroupTopology &gtopo = gc.GetGroupTopology();
3264  ReduceBeginCopy(x); // copy to 'ext_buf'
3265  int req_counter = 0;
3266  for (int nbr = 1; nbr < gtopo.GetNumNeighbors(); nbr++)
3267  {
3268  const int send_offset = ext_buf_offsets[nbr];
3269  const int send_size = ext_buf_offsets[nbr+1] - send_offset;
3270  if (send_size > 0)
3271  {
3272  auto send_buf = mpi_gpu_aware ? ext_buf.Read() : ext_buf.HostRead();
3273  MPI_Isend(send_buf + send_offset, send_size, MPI_DOUBLE,
3274  gtopo.GetNeighborRank(nbr), 41823,
3275  gtopo.GetComm(), &requests[req_counter++]);
3276  }
3277  const int recv_offset = shr_buf_offsets[nbr];
3278  const int recv_size = shr_buf_offsets[nbr+1] - recv_offset;
3279  if (recv_size > 0)
3280  {
3281  auto recv_buf = mpi_gpu_aware ? shr_buf.Write() : shr_buf.HostWrite();
3282  MPI_Irecv(recv_buf + recv_offset, recv_size, MPI_DOUBLE,
3283  gtopo.GetNeighborRank(nbr), 41823,
3284  gtopo.GetComm(), &requests[req_counter++]);
3285  }
3286  }
3287  ReduceLocalCopy(x, y);
3288  MPI_Waitall(req_counter, requests, MPI_STATUSES_IGNORE);
3289  ReduceEndAssemble(y); // assemble from 'shr_buf'
3290 }
3291 
3292 } // namespace mfem
3293 
3294 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
Geometry::Type GetGeometryType() const
Definition: element.hpp:52
int GetNFaceNeighbors() const
Definition: pmesh.hpp:286
int GetGroupMasterRank(int g) const
Return the rank of the group master for group &#39;g&#39;.
void SendRebalanceDofs(int old_ndofs, const Table &old_element_dofs, long old_global_offset, FiniteElementSpace *space)
Use the communication pattern from last Rebalance() to send element DOFs.
Definition: pncmesh.cpp:2133
void ReduceEnd(T *ldata, int layout, void(*Op)(OpData< T >)) const
Finalize reduction operation started with ReduceBegin().
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
void GetNeighborLDofTable(Table &nbr_ldof) const
Dofs to be received from communication neighbors.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
HYPRE_Int GetGlobalTDofNumber(int ldof) const
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:825
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:397
int * GetJ()
Definition: table.hpp:114
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
Definition: fespace.hpp:108
L2FaceValues
Definition: restriction.hpp:26
std::vector< int > CommGroup
Definition: pncmesh.hpp:131
Operator that extracts Face degrees of freedom in parallel.
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:144
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
Definition: ncmesh.cpp:4559
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T...
Definition: array.hpp:242
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
const CommGroup & GetGroup(GroupId id) const
Return a list of ranks contained in the group of the given ID.
Definition: pncmesh.hpp:160
void GetNeighborLTDofTable(Table &nbr_ltdof) const
Dofs to be sent to communication neighbors.
Operator that extracts Face degrees of freedom.
Array< int > ldof_group
Definition: nurbs.hpp:420
void BuildElementToDofTable() const
Definition: fespace.cpp:215
void AddColumnsInRow(int r, int ncol)
Definition: table.hpp:78
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:761
ConformingProlongationOperator(const ParFiniteElementSpace &pfes)
Definition: pfespace.cpp:2963
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
Definition: table.cpp:85
HYPRE_Int * GetDofOffsets() const
Definition: pfespace.hpp:247
static const int NumGeom
Definition: geom.hpp:41
int GetNGroups() const
Definition: pmesh.hpp:262
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
void ReduceEndAssemble(Vector &dst) const
Definition: pfespace.cpp:3252
Auxiliary class used by ParFiniteElementSpace.
Definition: pfespace.hpp:387
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
SparseMatrix * RefinementMatrix(int old_ndofs, const Table *old_elem_dof)
Definition: fespace.cpp:1095
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:173
Ordering::Type ordering
Definition: fespace.hpp:105
Array< Element * > face_nbr_elements
Definition: pmesh.hpp:254
void Delete()
Delete the owned pointers. The Memory is not reset by this method, i.e. it will, generally, not be Empty() after this call.
const Geometry::Type geom
Definition: ex1.cpp:40
void Create(const Array< int > &ldof_group)
Initialize the communicator from a local-dof to group map. Finalize() is called internally.
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:382
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:294
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:831
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
int GetGroupSize(int g) const
Get the number of processors in a group.
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2875
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
void SetDims(int rows, int nnz)
Definition: table.cpp:144
Operator::Type Type() const
Get the currently set operator type id.
Definition: handle.hpp:91
int VDofToDof(int vdof) const
Definition: fespace.hpp:507
void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:89
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
void ExchangeFaceNbrData(Table *gr_sface, int *s2l_face)
Definition: pmesh.cpp:1864
int GetNVertices() const
Definition: ncmesh.hpp:117
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs, int component=-1) const
Mark degrees of freedom associated with boundary elements with the specified boundary attributes (mar...
Definition: fespace.cpp:346
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:100
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:191
virtual void GetTrueTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer true-dof data from coarse_fes, defined on a coarse mesh, to this FE space, defined on a refined mesh.
Definition: pfespace.cpp:2848
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:715
void LoseData()
Call this if data has been stolen.
Definition: table.hpp:154
void BcastEndCopy(Vector &dst) const
Definition: pfespace.cpp:3164
void BcastBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3135
Abstract parallel finite element space.
Definition: pfespace.hpp:28
bool operator<(const Pair< A, B > &p, const Pair< A, B > &q)
Comparison operator for class Pair, based on the first element only.
Definition: sort_pairs.hpp:36
void Synchronize(Array< int > &ldof_marker) const
Given an integer array on the local degrees of freedom, perform a bitwise OR between the shared dofs...
Definition: pfespace.cpp:749
bool GroupContains(GroupId id, int rank) const
Return true if group &#39;id&#39; contains the given rank.
Definition: pncmesh.cpp:559
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:380
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:804
HypreParMatrix * GetPartialConformingInterpolation()
For a non-conforming mesh, construct and return the interpolation matrix from the partially conformin...
Definition: pfespace.cpp:704
const FiniteElementCollection * fec
Associated FE collection (not owned).
Definition: fespace.hpp:97
OperatorHandle Th
Transformation to apply to GridFunctions after space Update().
Definition: fespace.hpp:132
int GroupNQuadrilaterals(int group)
Definition: pmesh.hpp:268
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:834
int Size_of_connections() const
Definition: table.hpp:98
void GetLocalDerefinementMatrices(Geometry::Type geom, DenseTensor &localR) const
Definition: fespace.cpp:1379
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:885
bool HasGeometry(Geometry::Type geom) const
Return true iff the given geom is encountered in the mesh. Geometries of dimensions lower than Dimens...
Definition: mesh.hpp:846
void CopyRowStarts()
Definition: hypre.cpp:809
void DeleteAll()
Delete the whole array.
Definition: array.hpp:821
void AddConnections(int r, const int *c, int nc)
Definition: table.cpp:108
void GroupTriangle(int group, int i, int &face, int &o)
Definition: pmesh.cpp:1406
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
Definition: fespace.cpp:434
void BcastEnd(T *ldata, int layout) const
Finalize a broadcast started with BcastBegin().
bool IAmMaster(int g) const
Return true if I am master for group &#39;g&#39;.
void GetSharedEdgeDofs(int group, int ei, Array< int > &dofs) const
Definition: pfespace.cpp:531
ParNCMesh * pncmesh
Definition: pmesh.hpp:260
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Returns the indexes of the degrees of freedom for i&#39;th face including the dofs for the edges and the ...
Definition: fespace.cpp:1898
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:255
int GetNRanks() const
Definition: pmesh.hpp:237
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:474
int GroupVertex(int group, int i)
Definition: pmesh.hpp:270
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition: mesh.cpp:4184
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1165
ID for class SparseMatrix.
Definition: operator.hpp:240
DeviceConformingProlongationOperator(const ParFiniteElementSpace &pfes)
Definition: pfespace.cpp:3056
const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:376
long GetSequence() const
Definition: mesh.hpp:1227
FaceType
Definition: mesh.hpp:45
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1199
ID for the base class Operator, i.e. any type.
Definition: operator.hpp:239
Array< int > face_nbr_elements_offset
Definition: pmesh.hpp:252
void write(std::ostream &os, T value)
Write &#39;value&#39; to stream.
Definition: binaryio.hpp:30
static const int Dimension[NumGeom]
Definition: geom.hpp:46
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:726
GroupId GetEntityGroupId(int entity, int index)
Definition: pncmesh.hpp:148
int GetNumNeighbors() const
Return the number of neighbors including the local processor.
Memory< int > & GetJMemory()
Definition: table.hpp:119
void Clear()
Definition: table.cpp:381
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: pfespace.cpp:3260
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
MPI_Comm GetComm() const
Return the MPI communicator.
int GetNGhostEdges() const
Definition: pncmesh.hpp:103
void AddConnection(int r, int c)
Definition: table.hpp:80
void LoseData()
NULL-ifies the data.
Definition: array.hpp:118
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: fespace.hpp:330
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:1696
HYPRE_Int GetMyTDofOffset() const
Definition: pfespace.cpp:890
void ReduceLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3222
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
Definition: array.hpp:142
void SubDofOrder(Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:368
void GetSharedQuadrilateralDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:578
void CopyFrom(const U *src)
Copy from src into this array. Copies enough entries to fill the Capacity size of this array...
Definition: array.hpp:271
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:290
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:414
static void AddDependencies(SparseMatrix &deps, Array< int > &master_dofs, Array< int > &slave_dofs, DenseMatrix &I)
Definition: fespace.cpp:594
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition: fespace.hpp:710
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:951
T read(std::istream &is)
Read a value from the stream and return it.
Definition: binaryio.hpp:37
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:887
void Sort()
Sorts the array in ascending order. This requires operator&lt; to be defined for T.
Definition: array.hpp:234
int GroupNVertices(int group)
Definition: pmesh.hpp:265
void GetSharedTriangleDofs(int group, int fi, Array< int > &dofs) const
Definition: pfespace.cpp:554
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: pfespace.cpp:488
int FindSorted(const T &el) const
Do bisection search for &#39;el&#39; in a sorted array; return -1 if not found.
Definition: array.hpp:795
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:92
GroupTopology & GetGroupTopology()
Get a reference to the associated GroupTopology object.
int Dimension() const
Definition: mesh.hpp:788
int GetNFaces() const
Definition: ncmesh.hpp:119
Geometry::Type GetFaceGeometry(int index) const
Return face geometry type. index is the Mesh face number.
Definition: ncmesh.hpp:327
void GetEdgeDofs(int i, Array< int > &dofs) const
Returns the indexes of the degrees of freedom for i&#39;th edge including the dofs for the vertices of th...
Definition: fespace.cpp:1966
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
Auxiliary device class used by ParFiniteElementSpace.
Definition: pfespace.hpp:402
int GetGroupMaster(int g) const
Return the neighbor index of the group master for a given group. Neighbor 0 is the local processor...
void GetTransferOperator(const FiniteElementSpace &coarse_fes, OperatorHandle &T) const
Construct and return an Operator that can be used to transfer GridFunction data from coarse_fes...
Definition: fespace.cpp:2167
void GroupQuadrilateral(int group, int i, int &face, int &o)
Definition: pmesh.cpp:1417
int GroupNTriangles(int group)
Definition: pmesh.hpp:267
int GetNGhostFaces() const
Definition: pncmesh.hpp:104
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1031
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const =0
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:611
int GetMyRank() const
Definition: pmesh.hpp:238
Memory< int > & GetIMemory()
Definition: table.hpp:118
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1037
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: pfespace.cpp:3031
int GetNeighborRank(int i) const
Return the MPI rank of neighbor &#39;i&#39;.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void AddAColumnInRow(int r)
Definition: table.hpp:77
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:654
GroupId GetEntityOwnerId(int entity, int index)
Return vertex/edge/face (&#39;entity&#39; == 0/1/2, resp.) owner.
Definition: pncmesh.hpp:134
ParFiniteElementSpace(const ParFiniteElementSpace &orig, ParMesh *pmesh=NULL, const FiniteElementCollection *fec=NULL)
Copy constructor: deep copy all data from orig except the ParMesh, the FiniteElementCollection, and some derived data.
Definition: pfespace.cpp:30
HYPRE_Int GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:848
void SetLTDofTable(const Array< int > &ldof_ltdof)
Initialize the internal group_ltdof Table.
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition: device.hpp:246
Table & GroupLDofTable()
Fill-in the returned Table reference to initialize the GroupCommunicator then call Finalize()...
void GetEntityDofs(int entity, int index, Array< int > &dofs, Geometry::Type master_geom=Geometry::INVALID) const
Helper to get vertex, edge or face DOFs (entity=0,1,2 resp.).
Definition: fespace.cpp:670
int GetNEdges() const
Definition: ncmesh.hpp:118
void ShiftUpI()
Definition: table.cpp:119
Linear2DFiniteElement TriangleFE
Definition: fe.cpp:12833
virtual const Operator * GetFaceRestriction(ElementDofOrdering e_ordering, FaceType type, L2FaceValues mul=L2FaceValues::DoubleValued) const
Definition: pfespace.cpp:503
void BcastLocalCopy(const Vector &src, Vector &dst) const
Definition: pfespace.cpp:3156
void BcastBegin(T *ldata, int layout) const
Begin a broadcast within each group where the master is the root.
Operation GetLastOperation() const
Return type of last modification of the mesh.
Definition: mesh.hpp:1221
bool Nonconforming() const
Definition: pfespace.hpp:358
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1815
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:203
int GetMyRank() const
Definition: pncmesh.hpp:196
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: pfespace.cpp:460
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition: handle.hpp:112
Mesh * mesh
The mesh that FE space lives on (not owned).
Definition: fespace.hpp:94
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
GridFunction interpolation operator applicable after mesh refinement.
Definition: fespace.hpp:201
const NCList & GetNCList(int entity)
Return vertex/edge/face list (entity = 0/1/2, respectively).
Definition: ncmesh.hpp:240
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
Array< Master > masters
Definition: ncmesh.hpp:195
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
General triple product operator x -&gt; A*B*C*x, with ownership of the factors.
Definition: operator.hpp:769
void MakeJ()
Definition: table.cpp:95
int dim
Definition: ex24.cpp:53
Array< MeshId > conforming
Definition: ncmesh.hpp:194
int index(int i, int j, int nx, int ny)
Definition: life.cpp:241
bool IsGhost(int entity, int index) const
Return true if the specified vertex/edge/face is a ghost.
Definition: pncmesh.hpp:170
void CopyColStarts()
Definition: hypre.cpp:846
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3006
const FiniteElement * GetFaceNbrFaceFE(int i) const
Definition: pfespace.cpp:1213
double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:384
HYPRE_Int * GetTrueDofOffsets() const
Definition: pfespace.hpp:248
void ReduceBegin(const T *ldata) const
Begin reduction operation within each group where the master is the root.
void GetFaceNbrFaceVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1171
Geometry::Type GetFaceGeometryType(int Face) const
Definition: mesh.cpp:1043
void GetEdgeVertices(const MeshId &edge_id, int vert_index[2], bool oriented=true) const
Return Mesh vertex indices of an edge identified by &#39;edge_id&#39;.
Definition: ncmesh.cpp:4528
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:839
void GroupEdge(int group, int i, int &edge, int &o)
Definition: pmesh.cpp:1398
int GetNFaces() const
Return the number of faces in a 3D mesh.
Definition: mesh.hpp:746
BiLinear2DFiniteElement QuadrilateralFE
Vector data type.
Definition: vector.hpp:51
Identity Operator I: x -&gt; x.
Definition: operator.hpp:657
void GenerateOffsets(int N, HYPRE_Int loc_sizes[], Array< HYPRE_Int > *offsets[]) const
Definition: pmesh.cpp:1642
ID for class HypreParMatrix.
Definition: operator.hpp:241
int * GetI()
Definition: table.hpp:113
void ReduceBeginCopy(const Vector &src) const
Definition: pfespace.cpp:3211
static bool DofFinalizable(int dof, const Array< bool > &finalized, const SparseMatrix &deps)
Definition: fespace.cpp:618
int RowSize(int i) const
Definition: table.hpp:108
Biwise-OR of all device backends.
Definition: device.hpp:93
NURBSExtension * NURBSext
Definition: fespace.hpp:120
int GetNGhostVertices() const
Definition: pncmesh.hpp:102
Table send_face_nbr_elements
Definition: pmesh.hpp:257
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:118
GroupTopology gtopo
Definition: pmesh.hpp:247
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on scalar dofs, i.e. for VDim = 1.
Definition: pfespace.cpp:735
Class for parallel meshes.
Definition: pmesh.hpp:32
Abstract data type element.
Definition: element.hpp:28
GroupTopology gtopo
Definition: nurbs.hpp:418
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: pfespace.cpp:3172
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:2293
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:137
int GroupNEdges(int group)
Definition: pmesh.hpp:266
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.
double f(const Vector &p)
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:108
Array< HYPRE_Int > face_nbr_glob_dof_map
Definition: pfespace.hpp:184