MFEM  v3.1
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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "fem.hpp"
17 #include "../general/sort_pairs.hpp"
18 
19 #include <climits>
20 #include <list>
21 
22 namespace mfem
23 {
24 
25 ParFiniteElementSpace::ParFiniteElementSpace(ParFiniteElementSpace &pf)
26  : FiniteElementSpace(pf)
27 {
28  MyComm = pf.MyComm;
29  NRanks = pf.NRanks;
30  MyRank = pf.MyRank;
31  pmesh = pf.pmesh;
32  gcomm = pf.gcomm;
33  pf.gcomm = NULL;
34  ltdof_size = pf.ltdof_size;
35  Swap(ldof_group, pf.ldof_group);
36  Swap(ldof_ltdof, pf.ldof_ltdof);
37  Swap(dof_offsets, pf.dof_offsets);
38  Swap(tdof_offsets, pf.tdof_offsets);
39  Swap(tdof_nb_offsets, pf.tdof_nb_offsets);
40  Swap(ldof_sign, pf.ldof_sign);
41  P = pf.P;
42  pf.P = NULL;
43  R = pf.R;
44  pf.R = NULL;
45  num_face_nbr_dofs = pf.num_face_nbr_dofs;
46  pf.num_face_nbr_dofs = -1;
47  Swap<Table>(face_nbr_element_dof, pf.face_nbr_element_dof);
48  Swap<Table>(face_nbr_ldof, pf.face_nbr_ldof);
49  Swap(face_nbr_glob_dof_map, pf.face_nbr_glob_dof_map);
50  Swap<Table>(send_face_nbr_ldof, pf.send_face_nbr_ldof);
51 }
52 
53 ParFiniteElementSpace::ParFiniteElementSpace(
54  ParMesh *pm, const FiniteElementCollection *f, int dim, int ordering)
55  : FiniteElementSpace(pm, f, dim, ordering)
56 {
57  mesh = pmesh = pm;
58 
59  MyComm = pmesh->GetComm();
60  MPI_Comm_size(MyComm, &NRanks);
61  MPI_Comm_rank(MyComm, &MyRank);
62 
63  num_face_nbr_dofs = -1;
64 
65  P = NULL;
66  R = NULL;
67 
68  if (Nonconforming())
69  {
70  gcomm = NULL;
71  // needed here to initialize DOF offsets etc.
72  GetParallelConformingInterpolation();
73  return;
74  }
75 
76  if (NURBSext)
77  {
78  if (own_ext)
79  {
80  // the FiniteElementSpace constructor created a serial
81  // NURBSExtension of higher order than the mesh NURBSExtension
82 
84  NURBSext, dynamic_cast<ParNURBSExtension *>(pmesh->NURBSext));
85  // serial NURBSext is destroyed by the above constructor
86  NURBSext = pNe;
87  UpdateNURBS();
88  }
89 
90  ConstructTrueNURBSDofs();
91  }
92  else
93  {
94  ConstructTrueDofs();
95  }
96 
97  GenerateGlobalOffsets();
98 }
99 
100 void ParFiniteElementSpace::GetGroupComm(
101  GroupCommunicator &gc, int ldof_type, Array<int> *ldof_sign)
102 {
103  int gr;
104  int ng = pmesh->GetNGroups();
105  int nvd, ned, nfd;
106  Array<int> dofs;
107 
108  int group_ldof_counter;
109  Table &group_ldof = gc.GroupLDofTable();
110 
113  nfd = (fdofs) ? (fdofs[1]-fdofs[0]) : (0);
114 
115  if (ldof_sign)
116  {
117  ldof_sign->SetSize(GetNDofs());
118  *ldof_sign = 1;
119  }
120 
121  // count the number of ldofs in all groups (excluding the local group 0)
122  group_ldof_counter = 0;
123  for (gr = 1; gr < ng; gr++)
124  {
125  group_ldof_counter += nvd * pmesh->GroupNVertices(gr);
126  group_ldof_counter += ned * pmesh->GroupNEdges(gr);
127  group_ldof_counter += nfd * pmesh->GroupNFaces(gr);
128  }
129  if (ldof_type)
130  {
131  group_ldof_counter *= vdim;
132  }
133  // allocate the I and J arrays in group_ldof
134  group_ldof.SetDims(ng, group_ldof_counter);
135 
136  // build the full group_ldof table
137  group_ldof_counter = 0;
138  group_ldof.GetI()[0] = group_ldof.GetI()[1] = 0;
139  for (gr = 1; gr < ng; gr++)
140  {
141  int j, k, l, m, o, nv, ne, nf;
142  const int *ind;
143 
144  nv = pmesh->GroupNVertices(gr);
145  ne = pmesh->GroupNEdges(gr);
146  nf = pmesh->GroupNFaces(gr);
147 
148  // vertices
149  if (nvd > 0)
150  {
151  for (j = 0; j < nv; j++)
152  {
153  k = pmesh->GroupVertex(gr, j);
154 
155  dofs.SetSize(nvd);
156  m = nvd * k;
157  for (l = 0; l < nvd; l++, m++)
158  {
159  dofs[l] = m;
160  }
161 
162  if (ldof_type)
163  {
164  DofsToVDofs(dofs);
165  }
166 
167  for (l = 0; l < dofs.Size(); l++)
168  {
169  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
170  }
171  }
172  }
173 
174  // edges
175  if (ned > 0)
176  {
177  for (j = 0; j < ne; j++)
178  {
179  pmesh->GroupEdge(gr, j, k, o);
180 
181  dofs.SetSize(ned);
182  m = nvdofs+k*ned;
184  for (l = 0; l < ned; l++)
185  if (ind[l] < 0)
186  {
187  dofs[l] = m + (-1-ind[l]);
188  if (ldof_sign)
189  {
190  (*ldof_sign)[dofs[l]] = -1;
191  }
192  }
193  else
194  {
195  dofs[l] = m + ind[l];
196  }
197 
198  if (ldof_type)
199  {
200  DofsToVDofs(dofs);
201  }
202 
203  for (l = 0; l < dofs.Size(); l++)
204  {
205  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
206  }
207  }
208  }
209 
210  // faces
211  if (nfd > 0)
212  {
213  for (j = 0; j < nf; j++)
214  {
215  pmesh->GroupFace(gr, j, k, o);
216 
217  dofs.SetSize(nfd);
218  m = nvdofs+nedofs+fdofs[k];
220  mesh->GetFaceBaseGeometry(k), o);
221  for (l = 0; l < nfd; l++)
222  if (ind[l] < 0)
223  {
224  dofs[l] = m + (-1-ind[l]);
225  if (ldof_sign)
226  {
227  (*ldof_sign)[dofs[l]] = -1;
228  }
229  }
230  else
231  {
232  dofs[l] = m + ind[l];
233  }
234 
235  if (ldof_type)
236  {
237  DofsToVDofs(dofs);
238  }
239 
240  for (l = 0; l < dofs.Size(); l++)
241  {
242  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
243  }
244  }
245  }
246 
247  group_ldof.GetI()[gr+1] = group_ldof_counter;
248  }
249 
250  gc.Finalize();
251 }
252 
253 void ParFiniteElementSpace::ApplyLDofSigns(Array<int> &dofs) const
254 {
255  for (int i = 0; i < dofs.Size(); i++)
256  {
257  if (dofs[i] < 0)
258  {
259  if (ldof_sign[-1-dofs[i]] < 0)
260  {
261  dofs[i] = -1-dofs[i];
262  }
263  }
264  else
265  {
266  if (ldof_sign[dofs[i]] < 0)
267  {
268  dofs[i] = -1-dofs[i];
269  }
270  }
271  }
272 }
273 
275 {
276  if (elem_dof)
277  {
278  elem_dof->GetRow(i, dofs);
279  return;
280  }
282  if (Conforming())
283  {
284  ApplyLDofSigns(dofs);
285  }
286 }
287 
289 {
290  if (bdrElem_dof)
291  {
292  bdrElem_dof->GetRow(i, dofs);
293  return;
294  }
296  if (Conforming())
297  {
298  ApplyLDofSigns(dofs);
299  }
300 }
301 
303 {
305  if (Conforming())
306  {
307  ApplyLDofSigns(dofs);
308  }
309 }
310 
311 void ParFiniteElementSpace::GenerateGlobalOffsets()
312 {
313  if (HYPRE_AssumedPartitionCheck())
314  {
315  HYPRE_Int ldof[2];
316 
317  ldof[0] = GetVSize();
318  ldof[1] = TrueVSize();
319 
320  dof_offsets.SetSize(3);
321  tdof_offsets.SetSize(3);
322 
323  MPI_Scan(ldof, &dof_offsets[0], 2, HYPRE_MPI_INT, MPI_SUM, MyComm);
324 
325  tdof_offsets[1] = dof_offsets[1];
326  tdof_offsets[0] = tdof_offsets[1] - ldof[1];
327 
328  dof_offsets[1] = dof_offsets[0];
329  dof_offsets[0] = dof_offsets[1] - ldof[0];
330 
331  // get the global sizes in (t)dof_offsets[2]
332  if (MyRank == NRanks-1)
333  {
334  ldof[0] = dof_offsets[1];
335  ldof[1] = tdof_offsets[1];
336  }
337 
338  MPI_Bcast(ldof, 2, HYPRE_MPI_INT, NRanks-1, MyComm);
339  dof_offsets[2] = ldof[0];
340  tdof_offsets[2] = ldof[1];
341 
342  // Check for overflow
343  MFEM_VERIFY(dof_offsets[0] >= 0 && dof_offsets[1] >= 0,
344  "overflow in global dof_offsets");
345  MFEM_VERIFY(tdof_offsets[0] >= 0 && tdof_offsets[1] >= 0,
346  "overflow in global tdof_offsets");
347 
348  if (Conforming())
349  {
350  // Communicate the neighbor offsets in tdof_nb_offsets
351  GroupTopology &gt = GetGroupTopo();
352  int nsize = gt.GetNumNeighbors()-1;
353  MPI_Request *requests = new MPI_Request[2*nsize];
354  MPI_Status *statuses = new MPI_Status[2*nsize];
355  tdof_nb_offsets.SetSize(nsize+1);
356  tdof_nb_offsets[0] = tdof_offsets[0];
357 
358  int request_counter = 0;
359  // send and receive neighbors' local tdof offsets
360  for (int i = 1; i <= nsize; i++)
361  {
362  MPI_Irecv(&tdof_nb_offsets[i], 1, HYPRE_MPI_INT,
363  gt.GetNeighborRank(i), 5365, MyComm,
364  &requests[request_counter++]);
365  }
366  for (int i = 1; i <= nsize; i++)
367  {
368  MPI_Isend(&tdof_nb_offsets[0], 1, HYPRE_MPI_INT,
369  gt.GetNeighborRank(i), 5365, MyComm,
370  &requests[request_counter++]);
371  }
372  MPI_Waitall(request_counter, requests, statuses);
373 
374  delete [] statuses;
375  delete [] requests;
376  }
377  else
378  {
379  // TODO
380  }
381  }
382  else
383  {
384  int i;
385  HYPRE_Int ldof = GetVSize();
386  HYPRE_Int ltdof = TrueVSize();
387 
388  dof_offsets.SetSize (NRanks+1);
389  tdof_offsets.SetSize(NRanks+1);
390 
391  MPI_Allgather(&ldof, 1, HYPRE_MPI_INT,
392  &dof_offsets[1], 1, HYPRE_MPI_INT, MyComm);
393  MPI_Allgather(&ltdof, 1, HYPRE_MPI_INT,
394  &tdof_offsets[1], 1, HYPRE_MPI_INT, MyComm);
395 
396  dof_offsets[0] = tdof_offsets[0] = 0;
397  for (i = 1; i < NRanks; i++)
398  {
399  dof_offsets [i+1] += dof_offsets [i];
400  tdof_offsets[i+1] += tdof_offsets[i];
401  }
402 
403  // Check for overflow
404  MFEM_VERIFY(dof_offsets[MyRank] >= 0 && dof_offsets[MyRank+1] >= 0,
405  "overflow in global dof_offsets");
406  MFEM_VERIFY(tdof_offsets[MyRank] >= 0 && tdof_offsets[MyRank+1] >= 0,
407  "overflow in global tdof_offsets");
408  }
409 }
410 
412 {
413  if (P)
414  {
415  return P;
416  }
417 
418  if (Nonconforming())
419  {
420  GetParallelConformingInterpolation();
421  return P;
422  }
423 
424  int ldof = GetVSize();
425  int ltdof = TrueVSize();
426 
427  HYPRE_Int *i_diag = new HYPRE_Int[ldof+1];
428  HYPRE_Int *j_diag = new HYPRE_Int[ltdof];
429  int diag_counter;
430 
431  HYPRE_Int *i_offd = new HYPRE_Int[ldof+1];
432  HYPRE_Int *j_offd = new HYPRE_Int[ldof-ltdof];
433  int offd_counter;
434 
435  HYPRE_Int *cmap = new HYPRE_Int[ldof-ltdof];
436 
437  HYPRE_Int *col_starts = GetTrueDofOffsets();
438  HYPRE_Int *row_starts = GetDofOffsets();
439 
440  Array<Pair<HYPRE_Int, int> > cmap_j_offd(ldof-ltdof);
441 
442  i_diag[0] = i_offd[0] = 0;
443  diag_counter = offd_counter = 0;
444  for (int i = 0; i < ldof; i++)
445  {
446  int ltdof = GetLocalTDofNumber(i);
447  if (ltdof >= 0)
448  {
449  j_diag[diag_counter++] = ltdof;
450  }
451  else
452  {
453  cmap_j_offd[offd_counter].one = GetGlobalTDofNumber(i);
454  cmap_j_offd[offd_counter].two = offd_counter;
455  offd_counter++;
456  }
457  i_diag[i+1] = diag_counter;
458  i_offd[i+1] = offd_counter;
459  }
460 
461  SortPairs<HYPRE_Int, int>(cmap_j_offd, offd_counter);
462 
463  for (int i = 0; i < offd_counter; i++)
464  {
465  cmap[i] = cmap_j_offd[i].one;
466  j_offd[cmap_j_offd[i].two] = i;
467  }
468 
469  P = new HypreParMatrix(MyComm, MyRank, NRanks, row_starts, col_starts,
470  i_diag, j_diag, i_offd, j_offd, cmap, offd_counter);
471 
472  SparseMatrix Pdiag(0, 0);
473  P->GetDiag(Pdiag);
474  R = Transpose(Pdiag);
475 
476  return P;
477 }
478 
480 {
481  if (Nonconforming())
482  {
483  MFEM_ABORT("Not implemented for NC mesh.");
484  }
485 
486  GroupTopology &gt = GetGroupTopo();
487 
488  for (int i = 0; i < ldof_group.Size(); i++)
489  if (gt.IAmMaster(ldof_group[i])) // we are the master
490  {
491  vec[ldof_ltdof[i]] /= gt.GetGroupSize(ldof_group[i]);
492  }
493 }
494 
496 {
497  if (Nonconforming())
498  {
499  MFEM_WARNING("Not implemented for NC mesh.");
500  return NULL;
501  }
502 
503  GroupCommunicator *gc = new GroupCommunicator(GetGroupTopo());
504  if (NURBSext)
505  {
506  gc->Create(pNURBSext()->ldof_group);
507  }
508  else
509  {
510  GetGroupComm(*gc, 0);
511  }
512  return gc;
513 }
514 
516 {
517  if (Nonconforming())
518  {
519  MFEM_ABORT("Not implemented for NC mesh.");
520  }
521 
522  if (ldof_marker.Size() != GetVSize())
523  {
524  mfem_error("ParFiniteElementSpace::Synchronize");
525  }
526 
527  // implement allreduce(|) as reduce(|) + broadcast
528  gcomm->Reduce<int>(ldof_marker, GroupCommunicator::BitOR);
529  gcomm->Bcast(ldof_marker);
530 }
531 
533  Array<int> &ess_dofs) const
534 {
535  FiniteElementSpace::GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
536 
537  if (Conforming())
538  {
539  // Make sure that processors without boundary elements mark
540  // their boundary dofs (if they have any).
541  Synchronize(ess_dofs);
542  }
543 }
544 
546  const Array<int> &bdr_attr_is_ess, Array<int> &ess_tdof_list)
547 {
548  Array<int> ess_dofs, true_ess_dofs;
549 
550  GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
551  GetRestrictionMatrix()->BooleanMult(ess_dofs, true_ess_dofs);
552  MarkerToList(true_ess_dofs, ess_tdof_list);
553 }
554 
556 {
557  if (Nonconforming())
558  {
559  MFEM_VERIFY(P, "Dof_TrueDof_Matrix() needs to be called before "
560  "GetLocalTDofNumber()");
561 
562  return ldof_ltdof[ldof]; // NOTE: contains -1 for slaves/DOFs we don't own
563  }
564  else
565  {
566  if (GetGroupTopo().IAmMaster(ldof_group[ldof]))
567  {
568  return ldof_ltdof[ldof];
569  }
570  else
571  {
572  return -1;
573  }
574  }
575 }
576 
578 {
579  if (Nonconforming())
580  {
581  MFEM_VERIFY(ldof_ltdof[ldof] >= 0, "ldof " << ldof << " not a true DOF.");
582 
583  return GetMyTDofOffset() + ldof_ltdof[ldof];
584  }
585  else
586  {
587  if (HYPRE_AssumedPartitionCheck())
588  {
589  return ldof_ltdof[ldof] +
590  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(ldof_group[ldof])];
591  }
592  else
593  {
594  return ldof_ltdof[ldof] +
595  tdof_offsets[GetGroupTopo().GetGroupMasterRank(ldof_group[ldof])];
596  }
597  }
598 }
599 
601 {
602  if (Nonconforming())
603  {
604  MFEM_ABORT("Not implemented for NC mesh.");
605  }
606 
607  if (HYPRE_AssumedPartitionCheck())
608  {
610  {
611  return ldof_ltdof[sldof] +
612  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
613  ldof_group[sldof])] / vdim;
614  }
615  else
616  {
617  return (ldof_ltdof[sldof*vdim] +
618  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
619  ldof_group[sldof*vdim])]) / vdim;
620  }
621  }
622 
624  {
625  return ldof_ltdof[sldof] +
626  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
627  ldof_group[sldof])] / vdim;
628  }
629  else
630  {
631  return (ldof_ltdof[sldof*vdim] +
632  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
633  ldof_group[sldof*vdim])]) / vdim;
634  }
635 }
636 
638 {
639  return HYPRE_AssumedPartitionCheck() ? dof_offsets[0] : dof_offsets[MyRank];
640 }
641 
643 {
644  return HYPRE_AssumedPartitionCheck()? tdof_offsets[0] : tdof_offsets[MyRank];
645 }
646 
648 {
649  if (num_face_nbr_dofs >= 0) { return; }
650 
651  pmesh->ExchangeFaceNbrData();
652 
653  int num_face_nbrs = pmesh->GetNFaceNeighbors();
654 
655  if (num_face_nbrs == 0)
656  {
657  num_face_nbr_dofs = 0;
658  return;
659  }
660 
661  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
662  MPI_Request *send_requests = requests;
663  MPI_Request *recv_requests = requests + num_face_nbrs;
664  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
665 
666  Array<int> ldofs;
667  Array<int> ldof_marker(GetVSize());
668  ldof_marker = -1;
669 
670  Table send_nbr_elem_dof;
671 
672  send_nbr_elem_dof.MakeI(pmesh->send_face_nbr_elements.Size_of_connections());
673  send_face_nbr_ldof.MakeI(num_face_nbrs);
674  face_nbr_ldof.MakeI(num_face_nbrs);
675  int *send_el_off = pmesh->send_face_nbr_elements.GetI();
676  int *recv_el_off = pmesh->face_nbr_elements_offset;
677  for (int fn = 0; fn < num_face_nbrs; fn++)
678  {
679  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
680  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
681 
682  for (int i = 0; i < num_my_elems; i++)
683  {
684  GetElementVDofs(my_elems[i], ldofs);
685  for (int j = 0; j < ldofs.Size(); j++)
686  if (ldof_marker[ldofs[j]] != fn)
687  {
688  ldof_marker[ldofs[j]] = fn;
690  }
691  send_nbr_elem_dof.AddColumnsInRow(send_el_off[fn] + i, ldofs.Size());
692  }
693 
694  int nbr_rank = pmesh->GetFaceNbrRank(fn);
695  int tag = 0;
696  MPI_Isend(&send_face_nbr_ldof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
697  MyComm, &send_requests[fn]);
698 
699  MPI_Irecv(&face_nbr_ldof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
700  MyComm, &recv_requests[fn]);
701  }
702 
703  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
705 
707 
708  MPI_Waitall(num_face_nbrs, send_requests, statuses);
710 
711  // send/receive the I arrays of send_nbr_elem_dof/face_nbr_element_dof,
712  // respectively (they contain the number of dofs for each face-neighbor
713  // element)
714  face_nbr_element_dof.MakeI(recv_el_off[num_face_nbrs]);
715 
716  int *send_I = send_nbr_elem_dof.GetI();
717  int *recv_I = face_nbr_element_dof.GetI();
718  for (int fn = 0; fn < num_face_nbrs; fn++)
719  {
720  int nbr_rank = pmesh->GetFaceNbrRank(fn);
721  int tag = 0;
722  MPI_Isend(send_I + send_el_off[fn], send_el_off[fn+1] - send_el_off[fn],
723  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
724 
725  MPI_Irecv(recv_I + recv_el_off[fn], recv_el_off[fn+1] - recv_el_off[fn],
726  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
727  }
728 
729  MPI_Waitall(num_face_nbrs, send_requests, statuses);
730  send_nbr_elem_dof.MakeJ();
731 
732  ldof_marker = -1;
733 
734  for (int fn = 0; fn < num_face_nbrs; fn++)
735  {
736  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
737  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
738 
739  for (int i = 0; i < num_my_elems; i++)
740  {
741  GetElementVDofs(my_elems[i], ldofs);
742  for (int j = 0; j < ldofs.Size(); j++)
743  {
744  if (ldof_marker[ldofs[j]] != fn)
745  {
746  ldof_marker[ldofs[j]] = fn;
747  send_face_nbr_ldof.AddConnection(fn, ldofs[j]);
748  }
749  }
750  send_nbr_elem_dof.AddConnections(
751  send_el_off[fn] + i, ldofs, ldofs.Size());
752  }
753  }
755  send_nbr_elem_dof.ShiftUpI();
756 
757  // convert the ldof indices in send_nbr_elem_dof
758  int *send_J = send_nbr_elem_dof.GetJ();
759  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
760  {
761  int num_ldofs = send_face_nbr_ldof.RowSize(fn);
762  int *ldofs = send_face_nbr_ldof.GetRow(fn);
763  int j_end = send_I[send_el_off[fn+1]];
764 
765  for (int i = 0; i < num_ldofs; i++)
766  {
767  ldof_marker[ldofs[i]] = i;
768  }
769 
770  for ( ; j < j_end; j++)
771  {
772  send_J[j] = ldof_marker[send_J[j]];
773  }
774  }
775 
776  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
778 
779  // send/receive the J arrays of send_nbr_elem_dof/face_nbr_element_dof,
780  // respectively (they contain the element dofs in enumeration local for
781  // the face-neighbor pair)
782  int *recv_J = face_nbr_element_dof.GetJ();
783  for (int fn = 0; fn < num_face_nbrs; fn++)
784  {
785  int nbr_rank = pmesh->GetFaceNbrRank(fn);
786  int tag = 0;
787 
788  MPI_Isend(send_J + send_I[send_el_off[fn]],
789  send_I[send_el_off[fn+1]] - send_I[send_el_off[fn]],
790  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
791 
792  MPI_Irecv(recv_J + recv_I[recv_el_off[fn]],
793  recv_I[recv_el_off[fn+1]] - recv_I[recv_el_off[fn]],
794  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
795  }
796 
797  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
798 
799  // shift the J array of face_nbr_element_dof
800  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
801  {
802  int shift = face_nbr_ldof.GetI()[fn];
803  int j_end = recv_I[recv_el_off[fn+1]];
804 
805  for ( ; j < j_end; j++)
806  {
807  recv_J[j] += shift;
808  }
809  }
810 
811  MPI_Waitall(num_face_nbrs, send_requests, statuses);
812 
813  // send/receive the J arrays of send_face_nbr_ldof/face_nbr_ldof,
814  // respectively
815  send_J = send_face_nbr_ldof.GetJ();
816  for (int fn = 0; fn < num_face_nbrs; fn++)
817  {
818  int nbr_rank = pmesh->GetFaceNbrRank(fn);
819  int tag = 0;
820 
821  MPI_Isend(send_face_nbr_ldof.GetRow(fn),
823  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
824 
825  MPI_Irecv(face_nbr_ldof.GetRow(fn),
827  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
828  }
829 
830  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
831  MPI_Waitall(num_face_nbrs, send_requests, statuses);
832 
833  // send my_dof_offset (i.e. my_ldof_offset) to face neighbors and receive
834  // their offset in dof_face_nbr_offsets, used to define face_nbr_glob_dof_map
836  Array<HYPRE_Int> dof_face_nbr_offsets(num_face_nbrs);
837  HYPRE_Int my_dof_offset = GetMyDofOffset();
838  for (int fn = 0; fn < num_face_nbrs; fn++)
839  {
840  int nbr_rank = pmesh->GetFaceNbrRank(fn);
841  int tag = 0;
842 
843  MPI_Isend(&my_dof_offset, 1, HYPRE_MPI_INT, nbr_rank, tag,
844  MyComm, &send_requests[fn]);
845 
846  MPI_Irecv(&dof_face_nbr_offsets[fn], 1, HYPRE_MPI_INT, nbr_rank, tag,
847  MyComm, &recv_requests[fn]);
848  }
849 
850  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
851 
852  // set the array face_nbr_glob_dof_map which holds the global ldof indices of
853  // the face-neighbor dofs
854  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
855  {
856  for (int j_end = face_nbr_ldof.GetI()[fn+1]; j < j_end; j++)
858  dof_face_nbr_offsets[fn] + face_nbr_ldof.GetJ()[j];
859  }
860 
861  MPI_Waitall(num_face_nbrs, send_requests, statuses);
862 
863  delete [] statuses;
864  delete [] requests;
865 }
866 
868  int i, Array<int> &vdofs) const
869 {
870  face_nbr_element_dof.GetRow(i, vdofs);
871 }
872 
874 {
875  const FiniteElement *FE =
877  pmesh->face_nbr_elements[i]->GetGeometryType());
878 
879  if (NURBSext)
880  {
881  mfem_error("ParFiniteElementSpace::GetFaceNbrFE"
882  " does not support NURBS!");
883  }
884  return FE;
885 }
886 
888 {
889  hypre_ParCSRMatrix *csrP = (hypre_ParCSRMatrix*)(*P);
890  hypre_ParCSRMatrixOwnsRowStarts(csrP) = 1;
891  hypre_ParCSRMatrixOwnsColStarts(csrP) = 1;
892  P -> StealData();
893  dof_offsets.LoseData();
894  tdof_offsets.LoseData();
895 }
896 
897 void ParFiniteElementSpace::ConstructTrueDofs()
898 {
899  int i, gr, n = GetVSize();
900  GroupTopology &gt = pmesh->gtopo;
901  gcomm = new GroupCommunicator(gt);
902  Table &group_ldof = gcomm->GroupLDofTable();
903 
904  GetGroupComm(*gcomm, 1, &ldof_sign);
905 
906  // Define ldof_group and mark ldof_ltdof with
907  // -1 for ldof that is ours
908  // -2 for ldof that is in a group with another master
909  ldof_group.SetSize(n);
910  ldof_ltdof.SetSize(n);
911  ldof_group = 0;
912  ldof_ltdof = -1;
913 
914  for (gr = 1; gr < group_ldof.Size(); gr++)
915  {
916  const int *ldofs = group_ldof.GetRow(gr);
917  const int nldofs = group_ldof.RowSize(gr);
918  for (i = 0; i < nldofs; i++)
919  {
920  ldof_group[ldofs[i]] = gr;
921  }
922 
923  if (!gt.IAmMaster(gr)) // we are not the master
924  {
925  for (i = 0; i < nldofs; i++)
926  {
927  ldof_ltdof[ldofs[i]] = -2;
928  }
929  }
930  }
931 
932  // count ltdof_size
933  ltdof_size = 0;
934  for (i = 0; i < n; i++)
935  {
936  if (ldof_ltdof[i] == -1)
937  {
938  ldof_ltdof[i] = ltdof_size++;
939  }
940  }
941 
942  // have the group masters broadcast their ltdofs to the rest of the group
943  gcomm->Bcast(ldof_ltdof);
944 }
945 
946 void ParFiniteElementSpace::ConstructTrueNURBSDofs()
947 {
948  int n = GetVSize();
949  GroupTopology &gt = pNURBSext()->gtopo;
950  gcomm = new GroupCommunicator(gt);
951 
952  // pNURBSext()->ldof_group is for scalar space!
953  if (vdim == 1)
954  {
955  ldof_group.MakeRef(pNURBSext()->ldof_group);
956  }
957  else
958  {
959  const int *scalar_ldof_group = pNURBSext()->ldof_group;
960  ldof_group.SetSize(n);
961  for (int i = 0; i < n; i++)
962  {
963  ldof_group[i] = scalar_ldof_group[VDofToDof(i)];
964  }
965  }
966 
967  gcomm->Create(ldof_group);
968 
969  // ldof_sign.SetSize(n);
970  // ldof_sign = 1;
971  ldof_sign.DeleteAll();
972 
973  ltdof_size = 0;
974  ldof_ltdof.SetSize(n);
975  for (int i = 0; i < n; i++)
976  {
977  if (gt.IAmMaster(ldof_group[i]))
978  {
979  ldof_ltdof[i] = ltdof_size;
980  ltdof_size++;
981  }
982  else
983  {
984  ldof_ltdof[i] = -2;
985  }
986  }
987 
988  // have the group masters broadcast their ltdofs to the rest of the group
989  gcomm->Bcast(ldof_ltdof);
990 }
991 
992 inline int decode_dof(int dof, double& sign)
993 {
994  return (dof >= 0) ? (sign = 1, dof) : (sign = -1, (-1 - dof));
995 }
996 
997 const int INVALID_DOF = INT_MAX;
998 
999 static void MaskSlaveDofs(Array<int> &slave_dofs, const DenseMatrix &pm,
1000  const FiniteElementCollection *fec)
1001 {
1002  // If a slave edge/face shares a vertex with its master, that vertex is
1003  // not really a slave. In serial this is easily taken care of with the
1004  // statement "if (mdof != sdof)" inside AddDependencies. In parallel this
1005  // doesn't work as 'mdof' and 'sdof' may be on different processors.
1006  // Here we detect these cases from the slave point matrix.
1007 
1008  int k, i;
1009  int nv = fec->DofForGeometry(Geometry::POINT);
1010  int ne = fec->DofForGeometry(Geometry::SEGMENT);
1011 
1012  if (pm.Width() == 2) // edge: exclude master endpoint vertices
1013  {
1014  if (nv)
1015  {
1016  for (i = 0; i < 2; i++)
1017  {
1018  if (pm(0, i) == 0.0 || pm(0, i) == 1.0)
1019  {
1020  for (k = 0; k < nv; k++) { slave_dofs[i*nv + k] = INVALID_DOF; }
1021  }
1022  }
1023  }
1024  }
1025  else // face: exclude master corner vertices + full edges if present
1026  {
1027  MFEM_ASSERT(pm.Width() == 4, "");
1028 
1029  bool corner[4];
1030  for (i = 0; i < 4; i++)
1031  {
1032  double x = pm(0,i), y = pm(1,i);
1033  corner[i] = ((x == 0.0 || x == 1.0) && (y == 0.0 || y == 1.0));
1034  }
1035  if (nv)
1036  {
1037  for (i = 0; i < 4; i++)
1038  {
1039  if (corner[i])
1040  {
1041  for (k = 0; k < nv; k++) { slave_dofs[i*nv + k] = INVALID_DOF; }
1042  }
1043  }
1044  }
1045  if (ne)
1046  {
1047  for (i = 0; i < 4; i++)
1048  {
1049  if (corner[i] && corner[(i+1) % 4])
1050  {
1051  for (k = 0; k < ne; k++)
1052  {
1053  slave_dofs[4*nv + i*ne + k] = INVALID_DOF;
1054  }
1055  }
1056  }
1057  }
1058  }
1059 }
1060 
1061 void ParFiniteElementSpace
1062 ::AddSlaveDependencies(DepList deps[], int master_rank,
1063  const Array<int> &master_dofs, int master_ndofs,
1064  const Array<int> &slave_dofs, DenseMatrix& I)
1065 {
1066  // make each slave DOF dependent on all master DOFs
1067  for (int i = 0; i < slave_dofs.Size(); i++)
1068  {
1069  double ss, ms;
1070  int sdof = decode_dof(slave_dofs[i], ss);
1071  if (sdof == INVALID_DOF) { continue; }
1072 
1073  for (int vd = 0; vd < vdim; vd++)
1074  {
1075  DepList &dl = deps[DofToVDof(sdof, vd, ndofs)];
1076  if (dl.type < 2) // slave dependencies override 1-to-1 dependencies
1077  {
1078  Array<Dependency> tmp_list; // TODO remove, precalculate list size
1079  for (int j = 0; j < master_dofs.Size(); j++)
1080  {
1081  double coef = I(i, j);
1082  if (std::abs(coef) > 1e-12)
1083  {
1084  int mdof = decode_dof(master_dofs[j], ms);
1085  int mvdof = DofToVDof(mdof, vd, master_ndofs);
1086  tmp_list.Append(Dependency(master_rank, mvdof, coef*ms*ss));
1087  }
1088  }
1089  dl.type = 2;
1090  tmp_list.Copy(dl.list);
1091  }
1092  }
1093  }
1094 }
1095 
1096 void ParFiniteElementSpace
1097 ::Add1To1Dependencies(DepList deps[], int owner_rank,
1098  const Array<int> &owner_dofs, int owner_ndofs,
1099  const Array<int> &dependent_dofs)
1100 {
1101  MFEM_ASSERT(owner_dofs.Size() == dependent_dofs.Size(), "");
1102 
1103  for (int vd = 0; vd < vdim; vd++)
1104  {
1105  for (int i = 0; i < owner_dofs.Size(); i++)
1106  {
1107  double osign, dsign;
1108  int odof = decode_dof(owner_dofs[i], osign);
1109  int ddof = decode_dof(dependent_dofs[i], dsign);
1110  if (odof == INVALID_DOF || ddof == INVALID_DOF) { continue; }
1111 
1112  int ovdof = DofToVDof(odof, vd, owner_ndofs);
1113  int dvdof = DofToVDof(ddof, vd, ndofs);
1114 
1115  DepList &dl = deps[dvdof];
1116  if (dl.type == 0)
1117  {
1118  dl.type = 1;
1119  dl.list.Append(Dependency(owner_rank, ovdof, osign*dsign));
1120  }
1121  else if (dl.type == 1 && dl.list[0].rank > owner_rank)
1122  {
1123  // 1-to-1 dependency already exists but lower rank takes precedence
1124  dl.list[0] = Dependency(owner_rank, ovdof, osign*dsign);
1125  }
1126  }
1127  }
1128 }
1129 
1130 void ParFiniteElementSpace
1131 ::ReorderFaceDofs(Array<int> &dofs, int orient)
1132 {
1133  Array<int> tmp;
1134  dofs.Copy(tmp);
1135 
1136  // NOTE: assuming quad faces here
1137  int nv = fec->DofForGeometry(Geometry::POINT);
1139  const int* ind = fec->DofOrderForOrientation(Geometry::SQUARE, orient);
1140 
1141  int ve_dofs = 4*(nv + ne);
1142  for (int i = 0; i < ve_dofs; i++)
1143  {
1144  dofs[i] = INVALID_DOF;
1145  }
1146 
1147  int f_dofs = dofs.Size() - ve_dofs;
1148  for (int i = 0; i < f_dofs; i++)
1149  {
1150  if (ind[i] >= 0)
1151  {
1152  dofs[ve_dofs + i] = tmp[ve_dofs + ind[i]];
1153  }
1154  else
1155  {
1156  dofs[ve_dofs + i] = -1 - tmp[ve_dofs + (-1 - ind[i])];
1157  }
1158  }
1159 }
1160 
1161 void ParFiniteElementSpace::GetDofs(int type, int index, Array<int>& dofs)
1162 {
1163  switch (type)
1164  {
1165  case 0: GetVertexDofs(index, dofs); break;
1166  case 1: GetEdgeDofs(index, dofs); break;
1167  default: GetFaceDofs(index, dofs);
1168  }
1169 }
1170 
1171 void ParFiniteElementSpace::GetParallelConformingInterpolation()
1172 {
1173  ParNCMesh* pncmesh = pmesh->pncmesh;
1174 
1175  // *** STEP 1: exchange shared vertex/edge/face DOFs with neighbors ***
1176 
1177  NeighborDofMessage::Map send_dofs, recv_dofs;
1178 
1179  // prepare neighbor DOF messages for shared vertices/edges/faces
1180  for (int type = 0; type < 3; type++)
1181  {
1182  const NCMesh::NCList &list = pncmesh->GetSharedList(type);
1183  Array<int> dofs;
1184 
1185  int cs = list.conforming.size(), ms = list.masters.size();
1186  for (int i = 0; i < cs+ms; i++)
1187  {
1188  // loop through all (shared) conforming+master vertices/edges/faces
1189  const NCMesh::MeshId& id =
1190  (i < cs) ? (const NCMesh::MeshId&) list.conforming[i]
1191  /* */ : (const NCMesh::MeshId&) list.masters[i-cs];
1192 
1193  int owner = pncmesh->GetOwner(type, id.index), gsize;
1194  if (owner == MyRank)
1195  {
1196  // we own a shared v/e/f, send its DOFs to others in group
1197  GetDofs(type, id.index, dofs);
1198  // TODO: send dofs only when dofs.Size() > 0 ???
1199  const int *group = pncmesh->GetGroup(type, id.index, gsize);
1200  for (int j = 0; j < gsize; j++)
1201  {
1202  if (group[j] != MyRank)
1203  {
1204  NeighborDofMessage &send_msg = send_dofs[group[j]];
1205  send_msg.Init(pncmesh, fec, ndofs);
1206  send_msg.AddDofs(type, id, dofs);
1207  }
1208  }
1209  }
1210  else
1211  {
1212  // we don't own this v/e/f and expect to receive DOFs for it
1213  recv_dofs[owner].Init(pncmesh, fec, 0);
1214  }
1215  }
1216  }
1217 
1218  // send/receive all DOF messages
1219  NeighborDofMessage::IsendAll(send_dofs, MyComm);
1220  NeighborDofMessage::RecvAll(recv_dofs, MyComm);
1221 
1222  // *** STEP 2: build dependency lists ***
1223 
1224  int num_dofs = ndofs * vdim;
1225  DepList* deps = new DepList[num_dofs]; // NOTE: 'deps' is over vdofs
1226 
1227  Array<int> master_dofs, slave_dofs;
1228  Array<int> owner_dofs, my_dofs;
1229 
1230  // loop through *all* master edges/faces, constrain their slaves
1231  for (int type = 1; type < 3; type++)
1232  {
1233  const NCMesh::NCList &list = (type > 1) ? pncmesh->GetFaceList()
1234  /* */ : pncmesh->GetEdgeList();
1235  if (!list.masters.size()) { continue; }
1236 
1237  IsoparametricTransformation T;
1238  if (type > 1) { T.SetFE(&QuadrilateralFE); }
1239  else { T.SetFE(&SegmentFE); }
1240 
1241  int geom = (type > 1) ? Geometry::SQUARE : Geometry::SEGMENT;
1242  const FiniteElement* fe = fec->FiniteElementForGeometry(geom);
1243  if (!fe) { continue; }
1244 
1245  DenseMatrix I(fe->GetDof());
1246 
1247  // process masters that we own or that affect our edges/faces
1248  for (unsigned mi = 0; mi < list.masters.size(); mi++)
1249  {
1250  const NCMesh::Master &mf = list.masters[mi];
1251  if (!pncmesh->RankInGroup(type, mf.index, MyRank)) { continue; }
1252 
1253  // get master DOFs
1254  int master_ndofs, master_rank = pncmesh->GetOwner(type, mf.index);
1255  if (master_rank == MyRank)
1256  {
1257  GetDofs(type, mf.index, master_dofs);
1258  master_ndofs = ndofs;
1259  }
1260  else
1261  {
1262  recv_dofs[master_rank].GetDofs(type, mf, master_dofs, master_ndofs);
1263  }
1264 
1265  if (!master_dofs.Size()) { continue; }
1266 
1267  // constrain slaves that exist in our mesh
1268  for (int si = mf.slaves_begin; si < mf.slaves_end; si++)
1269  {
1270  const NCMesh::Slave &sf = list.slaves[si];
1271  if (pncmesh->IsGhost(type, sf.index)) { continue; }
1272 
1273  GetDofs(type, sf.index, slave_dofs);
1274  if (!slave_dofs.Size()) { continue; }
1275 
1276  T.GetPointMat() = sf.point_matrix;
1277  fe->GetLocalInterpolation(T, I);
1278 
1279  // make each slave DOF dependent on all master DOFs
1280  MaskSlaveDofs(slave_dofs, sf.point_matrix, fec);
1281  AddSlaveDependencies(deps, master_rank, master_dofs, master_ndofs,
1282  slave_dofs, I);
1283  }
1284 
1285  // special case for master edges that we don't own but still exist in
1286  // our mesh: this is a conforming-like situation, create 1-to-1 deps
1287  if (master_rank != MyRank && !pncmesh->IsGhost(type, mf.index))
1288  {
1289  GetDofs(type, mf.index, my_dofs);
1290  Add1To1Dependencies(deps, master_rank, master_dofs, master_ndofs,
1291  my_dofs);
1292  }
1293  }
1294  }
1295 
1296  // add one-to-one dependencies between shared conforming vertices/edges/faces
1297  for (int type = 0; type < 3; type++)
1298  {
1299  const NCMesh::NCList &list = pncmesh->GetSharedList(type);
1300  for (unsigned i = 0; i < list.conforming.size(); i++)
1301  {
1302  const NCMesh::MeshId &id = list.conforming[i];
1303  GetDofs(type, id.index, my_dofs);
1304  // TODO: skip if my_dofs.Size() == 0
1305 
1306  int owner_ndofs, owner = pncmesh->GetOwner(type, id.index);
1307  if (owner != MyRank)
1308  {
1309  recv_dofs[owner].GetDofs(type, id, owner_dofs, owner_ndofs);
1310  if (type == 2)
1311  {
1312  int fo = pncmesh->GetFaceOrientation(id.index);
1313  ReorderFaceDofs(owner_dofs, fo);
1314  }
1315  Add1To1Dependencies(deps, owner, owner_dofs, owner_ndofs, my_dofs);
1316  }
1317  else
1318  {
1319  // we own this v/e/f, assert ownership of the DOFs
1320  Add1To1Dependencies(deps, owner, my_dofs, ndofs, my_dofs);
1321  }
1322  }
1323  }
1324 
1325  // *** STEP 3: request P matrix rows that we need from neighbors ***
1326 
1327  NeighborRowRequest::Map send_requests, recv_requests;
1328 
1329  // copy communication topology from the DOF messages
1330  NeighborDofMessage::Map::iterator it;
1331  for (it = send_dofs.begin(); it != send_dofs.end(); ++it)
1332  {
1333  recv_requests[it->first];
1334  }
1335  for (it = recv_dofs.begin(); it != recv_dofs.end(); ++it)
1336  {
1337  send_requests[it->first];
1338  }
1339 
1340  // request rows we depend on
1341  for (int i = 0; i < num_dofs; i++)
1342  {
1343  const DepList &dl = deps[i];
1344  for (int j = 0; j < dl.list.Size(); j++)
1345  {
1346  const Dependency &dep = dl.list[j];
1347  if (dep.rank != MyRank)
1348  {
1349  send_requests[dep.rank].RequestRow(dep.dof);
1350  }
1351  }
1352  }
1353 
1354  NeighborRowRequest::IsendAll(send_requests, MyComm);
1355  NeighborRowRequest::RecvAll(recv_requests, MyComm);
1356 
1357  // *** STEP 4: iteratively build the P matrix ***
1358 
1359  // DOFs that stayed independent or are ours are true DOFs
1360  ltdof_size = 0;
1361  for (int i = 0; i < num_dofs; i++)
1362  {
1363  if (deps[i].IsTrueDof(MyRank)) { ltdof_size++; }
1364  }
1365 
1366  GenerateGlobalOffsets();
1367 
1368  HYPRE_Int glob_true_dofs = tdof_offsets.Last();
1369  HYPRE_Int glob_cdofs = dof_offsets.Last();
1370 
1371  // create the local part (local rows) of the P matrix
1372 #ifdef HYPRE_BIGINT
1373  MFEM_VERIFY(glob_true_dofs >= 0 && glob_true_dofs < (1ll << 31),
1374  "64bit matrix size not supported yet in non-conforming P.")
1375 #else
1376  MFEM_VERIFY(glob_true_dofs >= 0,
1377  "overflow of non-conforming P matrix columns.")
1378 #endif
1379  SparseMatrix localP(num_dofs, glob_true_dofs); // FIXME bigint
1380 
1381  // initialize the R matrix (also parallel but block-diagonal)
1382  R = new SparseMatrix(ltdof_size, num_dofs);
1383 
1384  Array<bool> finalized(num_dofs);
1385  finalized = false;
1386 
1387  // put identity in P and R for true DOFs, set ldof_ltdof
1388  HYPRE_Int my_tdof_offset = GetMyTDofOffset();
1389  ldof_ltdof.SetSize(num_dofs);
1390  ldof_ltdof = -1;
1391  for (int i = 0, true_dof = 0; i < num_dofs; i++)
1392  {
1393  if (deps[i].IsTrueDof(MyRank))
1394  {
1395  localP.Add(i, my_tdof_offset + true_dof, 1.0);
1396  R->Add(true_dof, i, 1.0);
1397  finalized[i] = true;
1398  ldof_ltdof[i] = true_dof;
1399  true_dof++;
1400  }
1401  }
1402 
1403  Array<int> cols;
1404  Vector srow;
1405 
1406  NeighborRowReply::Map recv_replies;
1407  std::list<NeighborRowReply::Map> send_replies;
1408 
1409  int num_finalized = ltdof_size;
1410  while (1)
1411  {
1412  // finalize all rows that can currently be finalized
1413  bool done;
1414  do
1415  {
1416  done = true;
1417  for (int dof = 0, i; dof < num_dofs; dof++)
1418  {
1419  if (finalized[dof]) { continue; }
1420 
1421  // check that rows of all constraining DOFs are available
1422  const DepList &dl = deps[dof];
1423  for (i = 0; i < dl.list.Size(); i++)
1424  {
1425  const Dependency &dep = dl.list[i];
1426  if (dep.rank == MyRank)
1427  {
1428  if (!finalized[dep.dof]) { break; }
1429  }
1430  else if (!recv_replies[dep.rank].HaveRow(dep.dof))
1431  {
1432  break;
1433  }
1434  }
1435  if (i < dl.list.Size()) { continue; }
1436 
1437  // form a linear combination of rows that 'dof' depends on
1438  for (i = 0; i < dl.list.Size(); i++)
1439  {
1440  const Dependency &dep = dl.list[i];
1441  if (dep.rank == MyRank)
1442  {
1443  localP.GetRow(dep.dof, cols, srow);
1444  }
1445  else
1446  {
1447  recv_replies[dep.rank].GetRow(dep.dof, cols, srow);
1448  }
1449  srow *= dep.coef;
1450  localP.AddRow(dof, cols, srow);
1451  }
1452 
1453  finalized[dof] = true;
1454  num_finalized++;
1455  done = false;
1456  }
1457  }
1458  while (!done);
1459 
1460  // send rows that are requested by neighbors and are available
1461  send_replies.push_back(NeighborRowReply::Map());
1462 
1463  NeighborRowRequest::Map::iterator it;
1464  for (it = recv_requests.begin(); it != recv_requests.end(); ++it)
1465  {
1466  NeighborRowRequest &req = it->second;
1467  std::set<int>::iterator row;
1468  for (row = req.rows.begin(); row != req.rows.end(); )
1469  {
1470  if (finalized[*row])
1471  {
1472  localP.GetRow(*row, cols, srow);
1473  send_replies.back()[it->first].AddRow(*row, cols, srow);
1474  req.rows.erase(row++);
1475  }
1476  else { ++row; }
1477  }
1478  }
1479  NeighborRowReply::IsendAll(send_replies.back(), MyComm);
1480 
1481  // are we finished?
1482  if (num_finalized >= num_dofs) { break; }
1483 
1484  // wait for a reply from neighbors
1485  int rank, size;
1486  NeighborRowReply::Probe(rank, size, MyComm);
1487  recv_replies[rank].Recv(rank, size, MyComm);
1488 
1489  // there may be more, receive all replies available
1490  while (NeighborRowReply::IProbe(rank, size, MyComm))
1491  {
1492  recv_replies[rank].Recv(rank, size, MyComm);
1493  }
1494  }
1495 
1496  delete [] deps;
1497  localP.Finalize();
1498 
1499  // create the parallel matrix P
1500 #ifndef HYPRE_BIGINT
1501  P = new HypreParMatrix(MyComm, num_dofs, glob_cdofs, glob_true_dofs,
1502  localP.GetI(), localP.GetJ(), localP.GetData(),
1503  dof_offsets.GetData(), tdof_offsets.GetData());
1504 #else
1505  (void) glob_cdofs;
1506  MFEM_ABORT("HYPRE_BIGINT not supported yet.");
1507 #endif
1508 
1509  R->Finalize();
1510 
1511  // make sure we can discard all send buffers
1513  NeighborRowRequest::WaitAllSent(send_requests);
1514 
1515  for (std::list<NeighborRowReply::Map>::iterator
1516  it = send_replies.begin(); it != send_replies.end(); ++it)
1517  {
1519  }
1520 }
1521 
1523 {
1525 
1526  ldof_group.DeleteAll();
1527  ldof_ltdof.DeleteAll();
1528  dof_offsets.DeleteAll();
1529  tdof_offsets.DeleteAll();
1530  tdof_nb_offsets.DeleteAll();
1531  ldof_sign.DeleteAll();
1532  delete P;
1533  P = NULL;
1534  delete R;
1535  R = NULL;
1536  delete gcomm;
1537  gcomm = NULL;
1538  num_face_nbr_dofs = -1;
1540  face_nbr_ldof.Clear();
1543  if (Conforming())
1544  {
1545  ConstructTrueDofs();
1546  GenerateGlobalOffsets();
1547  }
1548  else
1549  {
1550  GetParallelConformingInterpolation();
1551  }
1552 }
1553 
1555 {
1556  ParFiniteElementSpace *cpfes = new ParFiniteElementSpace(*this);
1557  Constructor();
1558  if (Conforming())
1559  {
1560  ConstructTrueDofs();
1561  GenerateGlobalOffsets();
1562  }
1563  else
1564  {
1565  GetParallelConformingInterpolation();
1566  }
1567  return cpfes;
1568 }
1569 
1570 }
1571 
1572 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:44
int GetNFaceNeighbors() const
Definition: pmesh.hpp:119
int GetGroupMasterRank(int g) const
int Size() const
Logical size of the array.
Definition: array.hpp:109
int GetVSize() const
Definition: fespace.hpp:164
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:1790
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:162
int * GetJ()
Definition: table.hpp:108
int ndofs
Number of degrees of freedom. Number of unknowns are ndofs*vdim.
Definition: fespace.hpp:72
template void SortPairs< HYPRE_Int, int >(Pair< HYPRE_Int, int > *, int)
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:116
Array< int > ldof_group
Definition: nurbs.hpp:350
void AddColumnsInRow(int r, int ncol)
Definition: table.hpp:72
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
Definition: table.cpp:74
void Swap< Table >(Table &a, Table &b)
Specialization of the template function Swap&lt;&gt; for class Table.
Definition: table.hpp:157
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:161
Array< Element * > face_nbr_elements
Definition: pmesh.hpp:97
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_vdofs) const
Definition: fespace.cpp:432
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:573
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
int GetGroupSize(int g) const
void SetDims(int rows, int nnz)
Definition: table.cpp:132
int VDofToDof(int vdof) const
Definition: fespace.hpp:251
T * GetData()
Returns the data.
Definition: array.hpp:91
HYPRE_Int * GetDofOffsets()
Definition: pfespace.hpp:151
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:69
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:179
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:479
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Synchronize(Array< int > &ldof_marker) const
Definition: pfespace.cpp:515
virtual const SparseMatrix * GetRestrictionMatrix()
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:215
std::map< int, NeighborRowRequest > Map
Definition: pncmesh.hpp:388
const FiniteElementCollection * fec
Definition: fespace.hpp:79
int Size_of_connections() const
Definition: table.hpp:92
HYPRE_Int GetMyDofOffset() const
Definition: pfespace.cpp:637
const int INVALID_DOF
Definition: pfespace.cpp:997
void DeleteAll()
Delete whole array.
Definition: array.hpp:455
void AddConnections(int r, const int *c, int nc)
Definition: table.cpp:96
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:487
bool IAmMaster(int g) const
ParNCMesh * pncmesh
Definition: pmesh.hpp:103
HYPRE_Int * GetTrueDofOffsets()
Definition: pfespace.hpp:152
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1295
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:288
int GroupVertex(int group, int i)
Definition: pmesh.hpp:112
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:867
void ExchangeFaceNbrData()
Definition: pmesh.cpp:843
void GetVertexDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1383
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:873
Array< int > face_nbr_elements_offset
Definition: pmesh.hpp:95
int dim
Definition: ex3.cpp:48
Data type sparse matrix.
Definition: sparsemat.hpp:38
void Clear()
Definition: table.cpp:340
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
void AddConnection(int r, int c)
Definition: table.hpp:74
void LoseData()
NULL-ifies the data.
Definition: array.hpp:103
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:1100
HYPRE_Int GetMyTDofOffset() const
Definition: pfespace.cpp:642
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:385
int GetLocalTDofNumber(int ldof)
Definition: pfespace.cpp:555
static void IsendAll(MapT &rank_msg, MPI_Comm comm)
Helper to send all messages in a rank-to-message map container.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:828
int GroupNVertices(int group)
Definition: pmesh.hpp:108
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: pfespace.cpp:302
int decode_dof(int dof, double &sign)
Definition: pfespace.cpp:992
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:86
virtual void Finalize(int skip_zeros=1)
Definition: sparsemat.cpp:680
void GetEdgeDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:1354
int GetGroupMaster(int g) const
HypreParMatrix * Dof_TrueDof_Matrix()
The true dof-to-dof interpolation matrix.
Definition: pfespace.cpp:411
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:314
virtual void Update()
Definition: fespace.cpp:1511
MPI_Comm GetComm() const
Definition: pmesh.hpp:86
Abstract finite element space.
Definition: fespace.hpp:62
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list)
Definition: pfespace.cpp:545
void AddAColumnInRow(int r)
Definition: table.hpp:71
void mfem_error(const char *msg)
Definition: error.cpp:23
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:323
static void RecvAll(MapT &rank_msg, MPI_Comm comm)
Helper to receive all messages in a rank-to-message map container.
HYPRE_Int GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:600
void ShiftUpI()
Definition: table.cpp:107
bool Nonconforming() const
Definition: pfespace.hpp:230
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: fespace.cpp:1212
NURBSExtension * NURBSext
Definition: mesh.hpp:142
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:274
void GroupFace(int group, int i, int &face, int &o)
Definition: pmesh.cpp:700
Mesh * mesh
The mesh that FE space lives on.
Definition: fespace.hpp:66
void MakeJ()
Definition: table.cpp:84
int GroupNFaces(int group)
Definition: pmesh.hpp:110
void Reduce(T *ldata, void(*Op)(OpData< T >))
void Create(Array< int > &ldof_group)
T & Last()
Return the last element in the array.
Definition: array.hpp:403
std::map< int, NeighborDofMessage > Map
Definition: pncmesh.hpp:360
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:466
void GroupEdge(int group, int i, int &edge, int &o)
Definition: pmesh.cpp:692
static void WaitAllSent(MapT &rank_msg)
Helper to wait for all messages in a map container to be sent.
void DofsToVDofs(Array< int > &dofs) const
Definition: fespace.cpp:35
int GetFaceBaseGeometry(int i) const
Definition: mesh.cpp:4348
int * GetI()
Definition: table.hpp:107
std::map< int, NeighborRowReply > Map
Definition: pncmesh.hpp:409
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:532
int RowSize(int i) const
Definition: table.hpp:102
NURBSExtension * NURBSext
Definition: fespace.hpp:90
virtual int DofForGeometry(int GeomType) const =0
static bool IProbe(int &rank, int &size, MPI_Comm comm)
Table send_face_nbr_elements
Definition: pmesh.hpp:100
GroupTopology gtopo
Definition: pmesh.hpp:90
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
int GetNGroups()
Definition: pmesh.hpp:105
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on Dofs.
Definition: pfespace.cpp:495
BiLinear2DFiniteElement QuadrilateralFE
Class for parallel meshes.
Definition: pmesh.hpp:28
virtual FiniteElementSpace * SaveUpdate()
Return a copy of the current FE space and update.
Definition: pfespace.cpp:1554
GroupTopology gtopo
Definition: nurbs.hpp:348
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:1293
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
HYPRE_Int GetGlobalTDofNumber(int ldof)
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:577
static void Probe(int &rank, int &size, MPI_Comm comm)
int GroupNEdges(int group)
Definition: pmesh.hpp:109
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.
Array< HYPRE_Int > face_nbr_glob_dof_map
Definition: pfespace.hpp:136