MFEM  v3.0
 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.googlecode.com.
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 namespace mfem
20 {
21 
22 ParFiniteElementSpace::ParFiniteElementSpace(ParFiniteElementSpace &pf)
23  : FiniteElementSpace(pf)
24 {
25  MyComm = pf.MyComm;
26  NRanks = pf.NRanks;
27  MyRank = pf.MyRank;
28  pmesh = pf.pmesh;
29  gcomm = pf.gcomm;
30  pf.gcomm = NULL;
31  ltdof_size = pf.ltdof_size;
32  Swap(ldof_group, pf.ldof_group);
33  Swap(ldof_ltdof, pf.ldof_ltdof);
34  Swap(dof_offsets, pf.dof_offsets);
35  Swap(tdof_offsets, pf.tdof_offsets);
36  Swap(tdof_nb_offsets, pf.tdof_nb_offsets);
37  Swap(ldof_sign, pf.ldof_sign);
38  P = pf.P;
39  pf.P = NULL;
40  num_face_nbr_dofs = pf.num_face_nbr_dofs;
41  pf.num_face_nbr_dofs = -1;
42  Swap<Table>(face_nbr_element_dof, pf.face_nbr_element_dof);
43  Swap<Table>(face_nbr_gdof, pf.face_nbr_gdof);
44  Swap<Table>(send_face_nbr_ldof, pf.send_face_nbr_ldof);
45 }
46 
47 ParFiniteElementSpace::ParFiniteElementSpace(
48  ParMesh *pm, const FiniteElementCollection *f, int dim, int order)
49  : FiniteElementSpace(pm, f, dim, order)
50 {
51  mesh = pmesh = pm;
52 
53  MyComm = pmesh->GetComm();
54  MPI_Comm_size(MyComm, &NRanks);
55  MPI_Comm_rank(MyComm, &MyRank);
56 
57  P = NULL;
58 
59  if (NURBSext)
60  {
61  if (own_ext)
62  {
63  // the FiniteElementSpace constructor created a serial
64  // NURBSExtension of higher order than the mesh NURBSExtension
65 
67  NURBSext, dynamic_cast<ParNURBSExtension *>(pmesh->NURBSext));
68  // serial NURBSext is destroyed by the above constructor
69  NURBSext = pNe;
70  UpdateNURBS();
71  }
72 
73  ConstructTrueNURBSDofs();
74  }
75  else
76  {
77  ConstructTrueDofs();
78  }
79 
80  GenerateGlobalOffsets();
81 
82  num_face_nbr_dofs = -1;
83 }
84 
85 void ParFiniteElementSpace::GetGroupComm(
86  GroupCommunicator &gc, int ldof_type, Array<int> *ldof_sign)
87 {
88  int gr;
89  int ng = pmesh->GetNGroups();
90  int nvd, ned, nfd;
91  Array<int> dofs;
92 
93  int group_ldof_counter;
94  Table &group_ldof = gc.GroupLDofTable();
95 
98  nfd = (fdofs) ? (fdofs[1]-fdofs[0]) : (0);
99 
100  if (ldof_sign)
101  {
102  ldof_sign->SetSize(GetNDofs());
103  *ldof_sign = 1;
104  }
105 
106  // count the number of ldofs in all groups (excluding the local group 0)
107  group_ldof_counter = 0;
108  for (gr = 1; gr < ng; gr++)
109  {
110  group_ldof_counter += nvd * pmesh->GroupNVertices(gr);
111  group_ldof_counter += ned * pmesh->GroupNEdges(gr);
112  group_ldof_counter += nfd * pmesh->GroupNFaces(gr);
113  }
114  if (ldof_type)
115  group_ldof_counter *= vdim;
116  // allocate the I and J arrays in group_ldof
117  group_ldof.SetDims(ng, group_ldof_counter);
118 
119  // build the full group_ldof table
120  group_ldof_counter = 0;
121  group_ldof.GetI()[0] = group_ldof.GetI()[1] = 0;
122  for (gr = 1; gr < ng; gr++)
123  {
124  int j, k, l, m, o, nv, ne, nf;
125  const int *ind;
126 
127  nv = pmesh->GroupNVertices(gr);
128  ne = pmesh->GroupNEdges(gr);
129  nf = pmesh->GroupNFaces(gr);
130 
131  // vertices
132  if (nvd > 0)
133  for (j = 0; j < nv; j++)
134  {
135  k = pmesh->GroupVertex(gr, j);
136 
137  dofs.SetSize(nvd);
138  m = nvd * k;
139  for (l = 0; l < nvd; l++, m++)
140  dofs[l] = m;
141 
142  if (ldof_type)
143  DofsToVDofs(dofs);
144 
145  for (l = 0; l < dofs.Size(); l++)
146  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
147  }
148 
149  // edges
150  if (ned > 0)
151  for (j = 0; j < ne; j++)
152  {
153  pmesh->GroupEdge(gr, j, k, o);
154 
155  dofs.SetSize(ned);
156  m = nvdofs+k*ned;
158  for (l = 0; l < ned; l++)
159  if (ind[l] < 0)
160  {
161  dofs[l] = m + (-1-ind[l]);
162  if (ldof_sign)
163  (*ldof_sign)[dofs[l]] = -1;
164  }
165  else
166  dofs[l] = m + ind[l];
167 
168  if (ldof_type)
169  DofsToVDofs(dofs);
170 
171  for (l = 0; l < dofs.Size(); l++)
172  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
173  }
174 
175  // faces
176  if (nfd > 0)
177  for (j = 0; j < nf; j++)
178  {
179  pmesh->GroupFace(gr, j, k, o);
180 
181  dofs.SetSize(nfd);
182  m = nvdofs+nedofs+fdofs[k];
184  mesh->GetFaceBaseGeometry(k), o);
185  for (l = 0; l < nfd; l++)
186  if (ind[l] < 0)
187  {
188  dofs[l] = m + (-1-ind[l]);
189  if (ldof_sign)
190  (*ldof_sign)[dofs[l]] = -1;
191  }
192  else
193  dofs[l] = m + ind[l];
194 
195  if (ldof_type)
196  DofsToVDofs(dofs);
197 
198  for (l = 0; l < dofs.Size(); l++)
199  group_ldof.GetJ()[group_ldof_counter++] = dofs[l];
200  }
201 
202  group_ldof.GetI()[gr+1] = group_ldof_counter;
203  }
204 
205  gc.Finalize();
206 }
207 
208 void ParFiniteElementSpace::ApplyLDofSigns(Array<int> &dofs) const
209 {
210  for (int i = 0; i < dofs.Size(); i++)
211  if (dofs[i] < 0)
212  {
213  if (ldof_sign[-1-dofs[i]] < 0)
214  dofs[i] = -1-dofs[i];
215  }
216  else
217  {
218  if (ldof_sign[dofs[i]] < 0)
219  dofs[i] = -1-dofs[i];
220  }
221 }
222 
224 {
225  if (elem_dof)
226  {
227  elem_dof->GetRow(i, dofs);
228  return;
229  }
231  ApplyLDofSigns(dofs);
232 }
233 
235 {
236  if (bdrElem_dof)
237  {
238  bdrElem_dof->GetRow(i, dofs);
239  return;
240  }
242  ApplyLDofSigns(dofs);
243 }
244 
246 {
248  ApplyLDofSigns(dofs);
249 }
250 
251 void ParFiniteElementSpace::GenerateGlobalOffsets()
252 {
253  if (HYPRE_AssumedPartitionCheck())
254  {
255  int ldof[2];
256 
257  ldof[0] = GetVSize();
258  ldof[1] = TrueVSize();
259 
260  dof_offsets.SetSize(3);
261  tdof_offsets.SetSize(3);
262 
263  MPI_Scan(ldof, &dof_offsets[0], 2, MPI_INT, MPI_SUM, MyComm);
264 
265  tdof_offsets[1] = dof_offsets[1];
266  tdof_offsets[0] = tdof_offsets[1] - ldof[1];
267 
268  dof_offsets[1] = dof_offsets[0];
269  dof_offsets[0] = dof_offsets[1] - ldof[0];
270 
271  // get the global sizes in (t)dof_offsets[2]
272  if (MyRank == NRanks-1)
273  {
274  ldof[0] = dof_offsets[1];
275  ldof[1] = tdof_offsets[1];
276  }
277 
278  MPI_Bcast(ldof, 2, MPI_INT, NRanks-1, MyComm);
279  dof_offsets[2] = ldof[0];
280  tdof_offsets[2] = ldof[1];
281  }
282  else
283  {
284  int i;
285  int ldof = GetVSize();
286  int ltdof = TrueVSize();
287 
288  dof_offsets.SetSize (NRanks+1);
289  tdof_offsets.SetSize(NRanks+1);
290 
291  MPI_Allgather(&ldof, 1, MPI_INT, &dof_offsets[1], 1, MPI_INT, MyComm);
292  MPI_Allgather(&ltdof, 1, MPI_INT, &tdof_offsets[1], 1, MPI_INT, MyComm);
293 
294  dof_offsets[0] = tdof_offsets[0] = 0;
295  for (i = 1; i < NRanks; i++)
296  {
297  dof_offsets [i+1] += dof_offsets [i];
298  tdof_offsets[i+1] += tdof_offsets[i];
299  }
300  }
301 }
302 
304 {
305  int i;
306 
307  if (P)
308  return P;
309 
310  int ldof = GetVSize();
311  int ltdof = TrueVSize();
312 
313  GroupTopology &gt = GetGroupTopo();
314 
315  int *i_diag;
316  int *j_diag;
317  int diag_counter;
318 
319  int *i_offd;
320  int *j_offd;
321  int offd_counter;
322 
323  int *cmap;
324  int *col_starts;
325  int *row_starts;
326 
327  col_starts = GetTrueDofOffsets();
328  row_starts = GetDofOffsets();
329 
330  i_diag = hypre_TAlloc(HYPRE_Int, ldof+1);
331  j_diag = hypre_TAlloc(HYPRE_Int, ltdof);
332 
333  i_offd = hypre_TAlloc(HYPRE_Int, ldof+1);
334  j_offd = hypre_TAlloc(HYPRE_Int, ldof-ltdof);
335 
336  cmap = hypre_TAlloc(HYPRE_Int, ldof-ltdof);
337 
338  Array<Pair<int, int> > cmap_j_offd(ldof-ltdof);
339 
340  if (HYPRE_AssumedPartitionCheck())
341  {
342  int nsize = gt.GetNumNeighbors()-1;
343  MPI_Request *requests = new MPI_Request[2*nsize];
344  MPI_Status *statuses = new MPI_Status[2*nsize];
345  tdof_nb_offsets.SetSize(nsize+1);
346  tdof_nb_offsets[0] = col_starts[0];
347 
348  int request_counter = 0;
349  // send and receive neighbors' local tdof offsets
350  for (i = 1; i <= nsize; i++)
351  MPI_Irecv(&tdof_nb_offsets[i], 1, MPI_INT, gt.GetNeighborRank(i), 5365,
352  MyComm, &requests[request_counter++]);
353 
354  for (i = 1; i <= nsize; i++)
355  MPI_Isend(&tdof_nb_offsets[0], 1, MPI_INT, gt.GetNeighborRank(i), 5365,
356  MyComm, &requests[request_counter++]);
357 
358  MPI_Waitall(request_counter, requests, statuses);
359 
360  delete [] statuses;
361  delete [] requests;
362  }
363 
364  i_diag[0] = i_offd[0] = 0;
365  diag_counter = offd_counter = 0;
366  for (i = 0; i < ldof; i++)
367  {
368  int proc = gt.GetGroupMasterRank(ldof_group[i]);
369  if (proc == MyRank)
370  {
371  j_diag[diag_counter++] = ldof_ltdof[i];
372  }
373  else
374  {
375  if (HYPRE_AssumedPartitionCheck())
376  cmap_j_offd[offd_counter].one =
377  tdof_nb_offsets[gt.GetGroupMaster(ldof_group[i])] + ldof_ltdof[i];
378  else
379  cmap_j_offd[offd_counter].one = col_starts[proc] + ldof_ltdof[i];
380  cmap_j_offd[offd_counter].two = offd_counter;
381  offd_counter++;
382  }
383  i_diag[i+1] = diag_counter;
384  i_offd[i+1] = offd_counter;
385  }
386 
387  SortPairs<int, int>(cmap_j_offd, offd_counter);
388 
389  for (i = 0; i < offd_counter; i++)
390  {
391  cmap[i] = cmap_j_offd[i].one;
392  j_offd[cmap_j_offd[i].two] = i;
393  }
394 
395  P = new HypreParMatrix(MyComm, MyRank, NRanks, row_starts, col_starts,
396  i_diag, j_diag, i_offd, j_offd, cmap, offd_counter);
397 
398  return P;
399 }
400 
402 {
403  GroupTopology &gt = GetGroupTopo();
404 
405  for (int i = 0; i < ldof_group.Size(); i++)
406  if (gt.IAmMaster(ldof_group[i])) // we are the master
407  vec[ldof_ltdof[i]] /= gt.GetGroupSize(ldof_group[i]);
408 }
409 
411 {
412  GroupCommunicator *gc = new GroupCommunicator(GetGroupTopo());
413  if (NURBSext)
414  gc->Create(pNURBSext()->ldof_group);
415  else
416  GetGroupComm(*gc, 0);
417  return gc;
418 }
419 
421 {
422  if (ldof_marker.Size() != GetVSize())
423  mfem_error("ParFiniteElementSpace::Synchronize");
424 
425  // implement allreduce(|) as reduce(|) + broadcast
426  gcomm->Reduce<int>(ldof_marker, GroupCommunicator::BitOR);
427  gcomm->Bcast(ldof_marker);
428 }
429 
431  Array<int> &ess_dofs) const
432 {
433  FiniteElementSpace::GetEssentialVDofs(bdr_attr_is_ess, ess_dofs);
434 
435  // Make sure that processors without boundary elements mark
436  // their boundary dofs (if they have any).
437  Synchronize(ess_dofs);
438 }
439 
441 {
442  if (GetGroupTopo().IAmMaster(ldof_group[ldof]))
443  return ldof_ltdof[ldof];
444  else
445  return -1;
446 }
447 
449 {
450  if (HYPRE_AssumedPartitionCheck())
451  {
452  if (!P)
454  return ldof_ltdof[ldof] +
455  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(ldof_group[ldof])];
456  }
457 
458  return ldof_ltdof[ldof] +
459  tdof_offsets[GetGroupTopo().GetGroupMasterRank(ldof_group[ldof])];
460 }
461 
463 {
464  if (HYPRE_AssumedPartitionCheck())
465  {
466  if (!P)
469  return ldof_ltdof[sldof] +
470  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
471  ldof_group[sldof])] / vdim;
472  else
473  return (ldof_ltdof[sldof*vdim] +
474  tdof_nb_offsets[GetGroupTopo().GetGroupMaster(
475  ldof_group[sldof*vdim])]) / vdim;
476  }
477 
479  return ldof_ltdof[sldof] +
480  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
481  ldof_group[sldof])] / vdim;
482  else
483  return (ldof_ltdof[sldof*vdim] +
484  tdof_offsets[GetGroupTopo().GetGroupMasterRank(
485  ldof_group[sldof*vdim])]) / vdim;
486 }
487 
489 {
490  if (HYPRE_AssumedPartitionCheck())
491  return dof_offsets[0];
492  else
493  return dof_offsets[MyRank];
494 }
495 
497 {
498  if (num_face_nbr_dofs >= 0)
499  return;
500 
501  pmesh->ExchangeFaceNbrData();
502 
503  int num_face_nbrs = pmesh->GetNFaceNeighbors();
504 
505  if (num_face_nbrs == 0)
506  {
507  num_face_nbr_dofs = 0;
508  return;
509  }
510 
511  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
512  MPI_Request *send_requests = requests;
513  MPI_Request *recv_requests = requests + num_face_nbrs;
514  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
515 
516  Array<int> ldofs;
517  Array<int> ldof_marker(GetVSize());
518  ldof_marker = -1;
519 
520  Table send_nbr_elem_dof;
521 
522  send_nbr_elem_dof.MakeI(pmesh->send_face_nbr_elements.Size_of_connections());
523  send_face_nbr_ldof.MakeI(num_face_nbrs);
524  face_nbr_gdof.MakeI(num_face_nbrs);
525  int *send_el_off = pmesh->send_face_nbr_elements.GetI();
526  int *recv_el_off = pmesh->face_nbr_elements_offset;
527  for (int fn = 0; fn < num_face_nbrs; fn++)
528  {
529  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
530  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
531 
532  for (int i = 0; i < num_my_elems; i++)
533  {
534  GetElementVDofs(my_elems[i], ldofs);
535  for (int j = 0; j < ldofs.Size(); j++)
536  if (ldof_marker[ldofs[j]] != fn)
537  {
538  ldof_marker[ldofs[j]] = fn;
540  }
541  send_nbr_elem_dof.AddColumnsInRow(send_el_off[fn] + i, ldofs.Size());
542  }
543 
544  int nbr_rank = pmesh->GetFaceNbrRank(fn);
545  int tag = 0;
546  MPI_Isend(&send_face_nbr_ldof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
547  MyComm, &send_requests[fn]);
548 
549  MPI_Irecv(&face_nbr_gdof.GetI()[fn], 1, MPI_INT, nbr_rank, tag,
550  MyComm, &recv_requests[fn]);
551  }
552 
553  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
555 
557 
558  MPI_Waitall(num_face_nbrs, send_requests, statuses);
560 
561  // send/receive the I arrays of send_nbr_elem_dof/face_nbr_element_dof,
562  // respectively (they contain the number of dofs for each face-neighbor
563  // element)
564  face_nbr_element_dof.MakeI(recv_el_off[num_face_nbrs]);
565 
566  int *send_I = send_nbr_elem_dof.GetI();
567  int *recv_I = face_nbr_element_dof.GetI();
568  for (int fn = 0; fn < num_face_nbrs; fn++)
569  {
570  int nbr_rank = pmesh->GetFaceNbrRank(fn);
571  int tag = 0;
572  MPI_Isend(send_I + send_el_off[fn], send_el_off[fn+1] - send_el_off[fn],
573  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
574 
575  MPI_Irecv(recv_I + recv_el_off[fn], recv_el_off[fn+1] - recv_el_off[fn],
576  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
577  }
578 
579  MPI_Waitall(num_face_nbrs, send_requests, statuses);
580  send_nbr_elem_dof.MakeJ();
581 
582  ldof_marker = -1;
583 
584  for (int fn = 0; fn < num_face_nbrs; fn++)
585  {
586  int *my_elems = pmesh->send_face_nbr_elements.GetRow(fn);
587  int num_my_elems = pmesh->send_face_nbr_elements.RowSize(fn);
588 
589  for (int i = 0; i < num_my_elems; i++)
590  {
591  GetElementVDofs(my_elems[i], ldofs);
592  for (int j = 0; j < ldofs.Size(); j++)
593  if (ldof_marker[ldofs[j]] != fn)
594  {
595  ldof_marker[ldofs[j]] = fn;
596  send_face_nbr_ldof.AddConnection(fn, ldofs[j]);
597  }
598  send_nbr_elem_dof.AddConnections(
599  send_el_off[fn] + i, ldofs, ldofs.Size());
600  }
601  }
603  send_nbr_elem_dof.ShiftUpI();
604 
605  // convert the ldof indices in send_nbr_elem_dof
606  int *send_J = send_nbr_elem_dof.GetJ();
607  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
608  {
609  int num_ldofs = send_face_nbr_ldof.RowSize(fn);
610  int *ldofs = send_face_nbr_ldof.GetRow(fn);
611  int j_end = send_I[send_el_off[fn+1]];
612 
613  for (int i = 0; i < num_ldofs; i++)
614  ldof_marker[ldofs[i]] = i;
615 
616  for ( ; j < j_end; j++)
617  send_J[j] = ldof_marker[send_J[j]];
618  }
619 
620  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
622 
623  // send/receive the J arrays of send_nbr_elem_dof/face_nbr_element_dof,
624  // respectively (they contain the element dofs in enumeration local for
625  // the face-neighbor pair)
626  int *recv_J = face_nbr_element_dof.GetJ();
627  for (int fn = 0; fn < num_face_nbrs; fn++)
628  {
629  int nbr_rank = pmesh->GetFaceNbrRank(fn);
630  int tag = 0;
631 
632  MPI_Isend(send_J + send_I[send_el_off[fn]],
633  send_I[send_el_off[fn+1]] - send_I[send_el_off[fn]],
634  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
635 
636  MPI_Irecv(recv_J + recv_I[recv_el_off[fn]],
637  recv_I[recv_el_off[fn+1]] - recv_I[recv_el_off[fn]],
638  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
639  }
640 
641  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
642 
643  // shift the J array of face_nbr_element_dof
644  for (int fn = 0, j = 0; fn < num_face_nbrs; fn++)
645  {
646  int shift = face_nbr_gdof.GetI()[fn];
647  int j_end = recv_I[recv_el_off[fn+1]];
648 
649  for ( ; j < j_end; j++)
650  recv_J[j] += shift;
651  }
652 
653  MPI_Waitall(num_face_nbrs, send_requests, statuses);
654 
655  // send/receive the J arrays of send_face_nbr_ldof/face_nbr_gdof,
656  // respectively
657  send_J = send_face_nbr_ldof.GetJ();
658  // switch to global dof numbers
659  int my_dof_offset = GetMyDofOffset();
660  int tot_send_dofs = send_face_nbr_ldof.Size_of_connections();
661  for (int i = 0; i < tot_send_dofs; i++)
662  send_J[i] += my_dof_offset;
663  for (int fn = 0; fn < num_face_nbrs; fn++)
664  {
665  int nbr_rank = pmesh->GetFaceNbrRank(fn);
666  int tag = 0;
667 
668  MPI_Isend(send_face_nbr_ldof.GetRow(fn),
670  MPI_INT, nbr_rank, tag, MyComm, &send_requests[fn]);
671 
672  MPI_Irecv(face_nbr_gdof.GetRow(fn),
674  MPI_INT, nbr_rank, tag, MyComm, &recv_requests[fn]);
675  }
676 
677  MPI_Waitall(num_face_nbrs, send_requests, statuses);
678 
679  // switch back to local dof numbers
680  for (int i = 0; i < tot_send_dofs; i++)
681  send_J[i] -= my_dof_offset;
682 
683  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
684 
685  delete [] statuses;
686  delete [] requests;
687 }
688 
690  int i, Array<int> &vdofs) const
691 {
692  face_nbr_element_dof.GetRow(i, vdofs);
693 }
694 
696 {
697  const FiniteElement *FE =
699  pmesh->face_nbr_elements[i]->GetGeometryType());
700 
701  if (NURBSext)
702  mfem_error("ParFiniteElementSpace::GetFaceNbrFE"
703  " does not support NURBS!");
704 
705  return FE;
706 }
707 
709 {
710  hypre_ParCSRMatrix *csrP = (hypre_ParCSRMatrix*)(*P);
711  hypre_ParCSRMatrixOwnsRowStarts(csrP) = 1;
712  hypre_ParCSRMatrixOwnsColStarts(csrP) = 1;
713  P -> StealData();
714  dof_offsets.LoseData();
715  tdof_offsets.LoseData();
716 }
717 
718 void ParFiniteElementSpace::ConstructTrueDofs()
719 {
720  int i, gr, n = GetVSize();
721  GroupTopology &gt = pmesh->gtopo;
722  gcomm = new GroupCommunicator(gt);
723  Table &group_ldof = gcomm->GroupLDofTable();
724 
725  GetGroupComm(*gcomm, 1, &ldof_sign);
726 
727  // Define ldof_group and mark ldof_ltdof with
728  // -1 for ldof that is ours
729  // -2 for ldof that is in a group with another master
730  ldof_group.SetSize(n);
731  ldof_ltdof.SetSize(n);
732  ldof_group = 0;
733  ldof_ltdof = -1;
734 
735  for (gr = 1; gr < group_ldof.Size(); gr++)
736  {
737  const int *ldofs = group_ldof.GetRow(gr);
738  const int nldofs = group_ldof.RowSize(gr);
739  for (i = 0; i < nldofs; i++)
740  ldof_group[ldofs[i]] = gr;
741 
742  if (!gt.IAmMaster(gr)) // we are not the master
743  for (i = 0; i < nldofs; i++)
744  ldof_ltdof[ldofs[i]] = -2;
745  }
746 
747  // count ltdof_size
748  ltdof_size = 0;
749  for (i = 0; i < n; i++)
750  if (ldof_ltdof[i] == -1)
751  ldof_ltdof[i] = ltdof_size++;
752 
753  // have the group masters broadcast their ltdofs to the rest of the group
754  gcomm->Bcast(ldof_ltdof);
755 }
756 
757 void ParFiniteElementSpace::ConstructTrueNURBSDofs()
758 {
759  int n = GetVSize();
760  GroupTopology &gt = pNURBSext()->gtopo;
761  gcomm = new GroupCommunicator(gt);
762 
763  // pNURBSext()->ldof_group is for scalar space!
764  if (vdim == 1)
765  {
766  ldof_group.MakeRef(pNURBSext()->ldof_group);
767  }
768  else
769  {
770  const int *scalar_ldof_group = pNURBSext()->ldof_group;
771  ldof_group.SetSize(n);
772  for (int i = 0; i < n; i++)
773  ldof_group[i] = scalar_ldof_group[VDofToDof(i)];
774  }
775 
776  gcomm->Create(ldof_group);
777 
778  // ldof_sign.SetSize(n);
779  // ldof_sign = 1;
780  ldof_sign.DeleteAll();
781 
782  ltdof_size = 0;
783  ldof_ltdof.SetSize(n);
784  for (int i = 0; i < n; i++)
785  {
786  if (gt.IAmMaster(ldof_group[i]))
787  {
788  ldof_ltdof[i] = ltdof_size;
789  ltdof_size++;
790  }
791  else
792  {
793  ldof_ltdof[i] = -2;
794  }
795  }
796 
797  // have the group masters broadcast their ltdofs to the rest of the group
798  gcomm->Bcast(ldof_ltdof);
799 }
800 
802 {
804 
805  ldof_group.DeleteAll();
806  ldof_ltdof.DeleteAll();
807  dof_offsets.DeleteAll();
808  tdof_offsets.DeleteAll();
809  tdof_nb_offsets.DeleteAll();
810  ldof_sign.DeleteAll();
811  delete P;
812  P = NULL;
813  delete gcomm;
814  gcomm = NULL;
815  num_face_nbr_dofs = -1;
819  ConstructTrueDofs();
820  GenerateGlobalOffsets();
821 }
822 
824 {
825  ParFiniteElementSpace *cpfes = new ParFiniteElementSpace(*this);
826  Constructor();
827  ConstructTrueDofs();
828  GenerateGlobalOffsets();
829  return cpfes;
830 }
831 
832 }
833 
834 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:42
int GetNFaceNeighbors() const
Definition: pmesh.hpp:99
int GetGroupMasterRank(int g) const
int Size() const
Logical size of the array.
Definition: array.hpp:108
int GetVSize() const
Definition: fespace.hpp:151
int GetNDofs() const
Returns number of degrees of freedom.
Definition: fespace.hpp:149
int * GetJ()
Definition: table.hpp:88
Array< int > ldof_group
Definition: nurbs.hpp:350
void AddColumnsInRow(int r, int ncol)
Definition: table.hpp:52
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
Definition: table.cpp:57
Array< Element * > face_nbr_elements
Definition: pmesh.hpp:79
MPI_Comm GetComm()
Definition: pmesh.hpp:68
int GetGroupSize(int g) const
void SetDims(int rows, int nnz)
Definition: table.cpp:107
int VDofToDof(int vdof) const
Definition: fespace.hpp:232
int vdim
Vector dimension (number of unknowns per degree of freedom).
Definition: fespace.hpp:68
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition: table.cpp:148
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:401
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Synchronize(Array< int > &ldof_marker) const
Definition: pfespace.cpp:420
const FiniteElementCollection * fec
Definition: fespace.hpp:78
int Size_of_connections() const
Definition: table.hpp:72
void DeleteAll()
Delete whole array.
Definition: array.hpp:407
void AddConnections(int r, const int *c, int nc)
Definition: table.cpp:75
bool IAmMaster(int g) const
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: fespace.cpp:888
template void SortPairs< int, int >(Pair< int, int > *, int)
virtual void GetBdrElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom for i&#39;th boundary element.
Definition: pfespace.cpp:234
int GroupVertex(int group, int i)
Definition: pmesh.hpp:92
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:689
void ExchangeFaceNbrData()
Definition: pmesh.cpp:662
int GetGlobalTDofNumber(int ldof)
Returns the global tdof number of the given local degree of freedom.
Definition: pfespace.cpp:448
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:695
Array< int > face_nbr_elements_offset
Definition: pmesh.hpp:77
int GetNumNeighbors() const
void Clear()
Definition: table.cpp:265
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
void AddConnection(int r, int c)
Definition: table.hpp:54
void LoseData()
NULL-ifies the data.
Definition: array.hpp:102
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:739
int GetLocalTDofNumber(int ldof)
Definition: pfespace.cpp:440
int GetGlobalScalarTDofNumber(int sldof)
Definition: pfespace.cpp:462
int GroupNVertices(int group)
Definition: pmesh.hpp:88
virtual void GetFaceDofs(int i, Array< int > &dofs) const
Definition: pfespace.cpp:245
int Size() const
Returns the number of TYPE I elements.
Definition: table.hpp:66
int GetGroupMaster(int g) const
HypreParMatrix * Dof_TrueDof_Matrix()
The true dof-to-dof interpolation matrix.
Definition: pfespace.cpp:303
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:284
virtual void Update()
Definition: fespace.cpp:1060
Abstract finite element space.
Definition: fespace.hpp:61
int GetNeighborRank(int i) const
void AddAColumnInRow(int r)
Definition: table.hpp:51
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:293
void GetElementVDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:117
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
Definition: fespace.cpp:362
void ShiftUpI()
Definition: table.cpp:84
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:825
NURBSExtension * NURBSext
Definition: mesh.hpp:307
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:223
void GroupFace(int group, int i, int &face, int &o)
Definition: pmesh.cpp:535
Mesh * mesh
The mesh that FE space lives on.
Definition: fespace.hpp:65
void MakeJ()
Definition: table.cpp:65
int GroupNFaces(int group)
Definition: pmesh.hpp:90
void Reduce(T *ldata, void(*Op)(OpData< T >))
void Create(Array< int > &ldof_group)
void MakeRef(T *, int)
Make this Array a reference to a poiter.
Definition: array.hpp:416
void GroupEdge(int group, int i, int &edge, int &o)
Definition: pmesh.cpp:527
void DofsToVDofs(Array< int > &dofs) const
Definition: fespace.cpp:29
int GetFaceBaseGeometry(int i) const
Definition: mesh.cpp:3589
int * GetI()
Definition: table.hpp:87
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs) const
Determine the boundary degrees of freedom.
Definition: pfespace.cpp:430
int RowSize(int i) const
Definition: table.hpp:82
NURBSExtension * NURBSext
Definition: fespace.hpp:89
virtual int DofForGeometry(int GeomType) const =0
Table send_face_nbr_elements
Definition: pmesh.hpp:82
GroupTopology gtopo
Definition: pmesh.hpp:72
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:103
int GetNGroups()
Definition: pmesh.hpp:85
GroupCommunicator * ScalarGroupComm()
Return a new GroupCommunicator on Dofs.
Definition: pfespace.cpp:410
Class for parallel meshes.
Definition: pmesh.hpp:27
virtual FiniteElementSpace * SaveUpdate()
Return a copy of the current FE space and update.
Definition: pfespace.cpp:823
GroupTopology gtopo
Definition: nurbs.hpp:348
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:1082
int GroupNEdges(int group)
Definition: pmesh.hpp:89
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
static void BitOR(OpData< T >)
Reduce operation bitwise OR, instantiated for int only.