MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pgridfunc.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 "fem.hpp"
17 #include <iostream>
18 #include <limits>
19 #include "../general/forall.hpp"
20 using namespace std;
21 
22 namespace mfem
23 {
24 
25 ParGridFunction::ParGridFunction(ParFiniteElementSpace *pf, GridFunction *gf)
26 {
27  fes = pfes = pf;
28  SetDataAndSize(gf->GetData(), gf->Size());
29 }
30 
31 ParGridFunction::ParGridFunction(ParFiniteElementSpace *pf, HypreParVector *tv)
32  : GridFunction(pf), pfes(pf)
33 {
34  Distribute(tv);
35 }
36 
38  const int *partitioning)
39 {
40  const FiniteElementSpace *glob_fes = gf->FESpace();
41  // duplicate the FiniteElementCollection from 'gf'
42  fec = FiniteElementCollection::New(glob_fes->FEColl()->Name());
43  // create a local ParFiniteElementSpace from the global one:
44  fes = pfes = new ParFiniteElementSpace(pmesh, glob_fes, partitioning, fec);
45  SetSize(pfes->GetVSize());
46 
47  if (partitioning)
48  {
49  // Assumption: the map "local element id" -> "global element id" is
50  // increasing, i.e. the local numbering preserves the element order from
51  // the global numbering.
52  Array<int> gvdofs, lvdofs;
53  Vector lnodes;
54  int element_counter = 0;
55  const int MyRank = pfes->GetMyRank();
56  const int glob_ne = glob_fes->GetNE();
57  for (int i = 0; i < glob_ne; i++)
58  {
59  if (partitioning[i] == MyRank)
60  {
61  pfes->GetElementVDofs(element_counter, lvdofs);
62  glob_fes->GetElementVDofs(i, gvdofs);
63  gf->GetSubVector(gvdofs, lnodes);
64  SetSubVector(lvdofs, lnodes);
65  element_counter++;
66  }
67  }
68  }
69 }
70 
71 ParGridFunction::ParGridFunction(ParMesh *pmesh, std::istream &input)
72  : GridFunction(pmesh, input)
73 {
74  // Convert the FiniteElementSpace, fes, to a ParFiniteElementSpace:
75  pfes = new ParFiniteElementSpace(pmesh, fec, fes->GetVDim(),
76  fes->GetOrdering());
77  delete fes;
78  fes = pfes;
79 }
80 
82 {
85 }
86 
88 {
91  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
92  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
93 }
94 
96 {
99  pfes = f;
100 }
101 
103 {
105  GridFunction::MakeRef(f, v);
106  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
107  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
108 }
109 
111 {
113  GridFunction::MakeRef(f, v);
114  pfes = f;
115 }
116 
118 {
120  GridFunction::MakeRef(f, v, v_offset);
121  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
122  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
123 }
124 
126 {
128  GridFunction::MakeRef(f, v, v_offset);
129  pfes = f;
130 }
131 
133 {
134  const Operator *prolong = pfes->GetProlongationMatrix();
135  prolong->Mult(*tv, *this);
136 }
137 
138 void ParGridFunction::AddDistribute(double a, const Vector *tv)
139 {
140  pfes->Dof_TrueDof_Matrix()->Mult(a, *tv, 1.0, *this);
141 }
142 
144 {
146  GetTrueDofs(*tv);
147  return tv;
148 }
149 
151 {
152  MFEM_VERIFY(pfes->Conforming(), "not implemented for NC meshes");
154  pfes->DivideByGroupSize(tv);
155 }
156 
158 {
159  MFEM_VERIFY(pfes->Conforming(), "not implemented for NC meshes");
161  pfes->DivideByGroupSize(tv);
162 }
163 
165 {
167  ParallelAverage(*tv);
168  return tv;
169 }
170 
172 {
173  pfes->GetRestrictionMatrix()->Mult(*this, tv);
174 }
175 
177 {
178  pfes->GetRestrictionMatrix()->Mult(*this, tv);
179 }
180 
182 {
184  ParallelProject(*tv);
185  return tv;
186 }
187 
189 {
191 }
192 
194 {
196 }
197 
199 {
201  ParallelAssemble(*tv);
202  return tv;
203 }
204 
206 {
208 
209  if (pfes->GetFaceNbrVSize() <= 0)
210  {
211  return;
212  }
213 
214  ParMesh *pmesh = pfes->GetParMesh();
215 
218 
219  int *send_offset = pfes->send_face_nbr_ldof.GetI();
220  const int *d_send_ldof = mfem::Read(pfes->send_face_nbr_ldof.GetJMemory(),
221  send_data.Size());
222  int *recv_offset = pfes->face_nbr_ldof.GetI();
223  MPI_Comm MyComm = pfes->GetComm();
224 
225  int num_face_nbrs = pmesh->GetNFaceNeighbors();
226  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
227  MPI_Request *send_requests = requests;
228  MPI_Request *recv_requests = requests + num_face_nbrs;
229  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
230 
231  auto d_data = this->Read();
232  auto d_send_data = send_data.Write();
233  MFEM_FORALL(i, send_data.Size(),
234  {
235  const int ldof = d_send_ldof[i];
236  d_send_data[i] = d_data[ldof >= 0 ? ldof : -1-ldof];
237  });
238 
239  bool mpi_gpu_aware = Device::GetGPUAwareMPI();
240  auto send_data_ptr = mpi_gpu_aware ? send_data.Read() : send_data.HostRead();
241  auto face_nbr_data_ptr = mpi_gpu_aware ? face_nbr_data.Write() :
243  for (int fn = 0; fn < num_face_nbrs; fn++)
244  {
245  int nbr_rank = pmesh->GetFaceNbrRank(fn);
246  int tag = 0;
247 
248  MPI_Isend(&send_data_ptr[send_offset[fn]],
249  send_offset[fn+1] - send_offset[fn],
250  MPI_DOUBLE, nbr_rank, tag, MyComm, &send_requests[fn]);
251 
252  MPI_Irecv(&face_nbr_data_ptr[recv_offset[fn]],
253  recv_offset[fn+1] - recv_offset[fn],
254  MPI_DOUBLE, nbr_rank, tag, MyComm, &recv_requests[fn]);
255  }
256 
257  MPI_Waitall(num_face_nbrs, send_requests, statuses);
258  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
259 
260  delete [] statuses;
261  delete [] requests;
262 }
263 
264 double ParGridFunction::GetValue(int i, const IntegrationPoint &ip, int vdim)
265 const
266 {
267  Array<int> dofs;
268  Vector DofVal, LocVec;
269  int nbr_el_no = i - pfes->GetParMesh()->GetNE();
270  if (nbr_el_no >= 0)
271  {
272  int fes_vdim = pfes->GetVDim();
273  pfes->GetFaceNbrElementVDofs(nbr_el_no, dofs);
274  const FiniteElement *fe = pfes->GetFaceNbrFE(nbr_el_no);
275  if (fes_vdim > 1)
276  {
277  int s = dofs.Size()/fes_vdim;
278  Array<int> _dofs(&dofs[(vdim-1)*s], s);
279  face_nbr_data.GetSubVector(_dofs, LocVec);
280  DofVal.SetSize(s);
281  }
282  else
283  {
284  face_nbr_data.GetSubVector(dofs, LocVec);
285  DofVal.SetSize(dofs.Size());
286  }
287  if (fe->GetMapType() == FiniteElement::VALUE)
288  {
289  fe->CalcShape(ip, DofVal);
290  }
291  else
292  {
295  Tr->SetIntPoint(&ip);
296  fe->CalcPhysShape(*Tr, DofVal);
297  }
298  }
299  else
300  {
301  fes->GetElementDofs(i, dofs);
302  fes->DofsToVDofs(vdim-1, dofs);
303  DofVal.SetSize(dofs.Size());
304  const FiniteElement *fe = fes->GetFE(i);
305  if (fe->GetMapType() == FiniteElement::VALUE)
306  {
307  fe->CalcShape(ip, DofVal);
308  }
309  else
310  {
312  Tr->SetIntPoint(&ip);
313  fe->CalcPhysShape(*Tr, DofVal);
314  }
315  GetSubVector(dofs, LocVec);
316  }
317 
318  return (DofVal * LocVec);
319 }
320 
322  Vector &val) const
323 {
324  int nbr_el_no = i - pfes->GetParMesh()->GetNE();
325  if (nbr_el_no >= 0)
326  {
327  Array<int> dofs;
328  pfes->GetFaceNbrElementVDofs(nbr_el_no, dofs);
329  Vector loc_data;
330  face_nbr_data.GetSubVector(dofs, loc_data);
331  const FiniteElement *FElem = pfes->GetFaceNbrFE(nbr_el_no);
332  int dof = FElem->GetDof();
333  if (FElem->GetRangeType() == FiniteElement::SCALAR)
334  {
335  Vector shape(dof);
336  if (FElem->GetMapType() == FiniteElement::VALUE)
337  {
338  FElem->CalcShape(ip, shape);
339  }
340  else
341  {
344  Tr->SetIntPoint(&ip);
345  FElem->CalcPhysShape(*Tr, shape);
346  }
347  int vdim = fes->GetVDim();
348  val.SetSize(vdim);
349  for (int k = 0; k < vdim; k++)
350  {
351  val(k) = shape * ((const double *)loc_data + dof * k);
352  }
353  }
354  else
355  {
356  int spaceDim = fes->GetMesh()->SpaceDimension();
357  DenseMatrix vshape(dof, spaceDim);
360  Tr->SetIntPoint(&ip);
361  FElem->CalcVShape(*Tr, vshape);
362  val.SetSize(spaceDim);
363  vshape.MultTranspose(loc_data, val);
364  }
365  }
366  else
367  {
368  GridFunction::GetVectorValue(i, ip, val);
369  }
370 }
371 
373  const IntegrationPoint &ip,
374  int comp, Vector *tr) const
375 {
376  // We can assume faces and edges are local
378  {
379  return GridFunction::GetValue(T, ip, comp, tr);
380  }
381 
382  // Check for evaluation in a local element
383  int nbr_el_no = T.ElementNo - pfes->GetParMesh()->GetNE();
384  if (nbr_el_no < 0)
385  {
386  return GridFunction::GetValue(T, ip, comp, tr);
387  }
388 
389  // Evaluate using DoFs from a neighboring element
390  if (tr)
391  {
392  T.SetIntPoint(&ip);
393  T.Transform(ip, *tr);
394  }
395 
396  Array<int> dofs;
397  const FiniteElement * fe = pfes->GetFaceNbrFE(nbr_el_no);
398  pfes->GetFaceNbrElementVDofs(nbr_el_no, dofs);
399 
400  pfes->DofsToVDofs(comp-1, dofs);
401  Vector DofVal(dofs.Size()), LocVec;
402  if (fe->GetMapType() == FiniteElement::VALUE)
403  {
404  fe->CalcShape(ip, DofVal);
405  }
406  else
407  {
408  fe->CalcPhysShape(T, DofVal);
409  }
410  face_nbr_data.GetSubVector(dofs, LocVec);
411 
412  return (DofVal * LocVec);
413 }
414 
416  const IntegrationPoint &ip,
417  Vector &val, Vector *tr) const
418 {
419  // We can assume faces and edges are local
421  {
422  return GridFunction::GetVectorValue(T, ip, val, tr);
423  }
424 
425  // Check for evaluation in a local element
426  int nbr_el_no = T.ElementNo - pfes->GetParMesh()->GetNE();
427  if (nbr_el_no < 0)
428  {
429  return GridFunction::GetVectorValue(T, ip, val, tr);
430  }
431 
432  // Evaluate using DoFs from a neighboring element
433  if (tr)
434  {
435  T.SetIntPoint(&ip);
436  T.Transform(ip, *tr);
437  }
438 
439  Array<int> vdofs;
440  pfes->GetFaceNbrElementVDofs(nbr_el_no, vdofs);
441  const FiniteElement *fe = pfes->GetFaceNbrFE(nbr_el_no);
442 
443  int dof = fe->GetDof();
444  Vector loc_data;
445  face_nbr_data.GetSubVector(vdofs, loc_data);
446  if (fe->GetRangeType() == FiniteElement::SCALAR)
447  {
448  Vector shape(dof);
449  if (fe->GetMapType() == FiniteElement::VALUE)
450  {
451  fe->CalcShape(ip, shape);
452  }
453  else
454  {
455  fe->CalcPhysShape(T, shape);
456  }
457  int vdim = pfes->GetVDim();
458  val.SetSize(vdim);
459  for (int k = 0; k < vdim; k++)
460  {
461  val(k) = shape * ((const double *)loc_data + dof * k);
462  }
463  }
464  else
465  {
466  int spaceDim = pfes->GetMesh()->SpaceDimension();
467  DenseMatrix vshape(dof, spaceDim);
468  fe->CalcVShape(T, vshape);
469  val.SetSize(spaceDim);
470  vshape.MultTranspose(loc_data, val);
471  }
472 }
473 
475 {
476  DeltaCoefficient *delta_c = dynamic_cast<DeltaCoefficient *>(&coeff);
477 
478  if (delta_c == NULL)
479  {
481  }
482  else
483  {
484  double loc_integral, glob_integral;
485 
486  ProjectDeltaCoefficient(*delta_c, loc_integral);
487 
488  MPI_Allreduce(&loc_integral, &glob_integral, 1, MPI_DOUBLE, MPI_SUM,
489  pfes->GetComm());
490 
491  (*this) *= (delta_c->Scale() / glob_integral);
492  }
493 }
494 
496 {
497  // local maximal element attribute for each dof
498  Array<int> ldof_attr;
499 
500  // local projection
501  GridFunction::ProjectDiscCoefficient(coeff, ldof_attr);
502 
503  // global maximal element attribute for each dof
504  Array<int> gdof_attr;
505  ldof_attr.Copy(gdof_attr);
506  GroupCommunicator &gcomm = pfes->GroupComm();
507  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Max);
508  gcomm.Bcast(gdof_attr);
509 
510  // set local value to zero if global maximal element attribute is larger than
511  // the local one, and mark (in gdof_attr) if we have the correct value
512  for (int i = 0; i < pfes->GetVSize(); i++)
513  {
514  if (gdof_attr[i] > ldof_attr[i])
515  {
516  (*this)(i) = 0.0;
517  gdof_attr[i] = 0;
518  }
519  else
520  {
521  gdof_attr[i] = 1;
522  }
523  }
524 
525  // parallel averaging plus interpolation to determine final values
527  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Sum);
528  gcomm.Bcast(gdof_attr);
529  for (int i = 0; i < fes->GetVSize(); i++)
530  {
531  (*this)(i) /= gdof_attr[i];
532  }
533  this->ParallelAssemble(*tv);
534  this->Distribute(tv);
535  delete tv;
536 }
537 
538 
540 {
541  // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
542  // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
543 
544  // Number of zones that contain a given dof.
545  Array<int> zones_per_vdof;
546  AccumulateAndCountZones(coeff, type, zones_per_vdof);
547 
548  // Count the zones globally.
549  GroupCommunicator &gcomm = pfes->GroupComm();
550  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
551  gcomm.Bcast(zones_per_vdof);
552 
553  // Accumulate for all vdofs.
554  gcomm.Reduce<double>(data, GroupCommunicator::Sum);
555  gcomm.Bcast<double>(data);
556 
557  ComputeMeans(type, zones_per_vdof);
558 }
559 
561  AvgType type)
562 {
563  // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
564  // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
565 
566  // Number of zones that contain a given dof.
567  Array<int> zones_per_vdof;
568  AccumulateAndCountZones(vcoeff, type, zones_per_vdof);
569 
570  // Count the zones globally.
571  GroupCommunicator &gcomm = pfes->GroupComm();
572  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
573  gcomm.Bcast(zones_per_vdof);
574 
575  // Accumulate for all vdofs.
576  gcomm.Reduce<double>(data, GroupCommunicator::Sum);
577  gcomm.Bcast<double>(data);
578 
579  ComputeMeans(type, zones_per_vdof);
580 }
581 
583  Coefficient *coeff[], VectorCoefficient *vcoeff, Array<int> &attr)
584 {
585  Array<int> values_counter;
586  AccumulateAndCountBdrValues(coeff, vcoeff, attr, values_counter);
587 
588  Vector values(Size());
589  for (int i = 0; i < values.Size(); i++)
590  {
591  values(i) = values_counter[i] ? (*this)(i) : 0.0;
592  }
593 
594  // Count the values globally.
595  GroupCommunicator &gcomm = pfes->GroupComm();
596  gcomm.Reduce<int>(values_counter, GroupCommunicator::Sum);
597  // Accumulate the values globally.
598  gcomm.Reduce<double>(values, GroupCommunicator::Sum);
599  // Only the values in the master are guaranteed to be correct!
600  for (int i = 0; i < values.Size(); i++)
601  {
602  if (values_counter[i])
603  {
604  (*this)(i) = values(i)/values_counter[i];
605  }
606  }
607 
608 #ifdef MFEM_DEBUG
609  Array<int> ess_vdofs_marker;
610  pfes->GetEssentialVDofs(attr, ess_vdofs_marker);
611  for (int i = 0; i < values_counter.Size(); i++)
612  {
613  MFEM_ASSERT(pfes->GetLocalTDofNumber(i) == -1 ||
614  bool(values_counter[i]) == bool(ess_vdofs_marker[i]),
615  "internal error");
616  }
617 #endif
618 }
619 
621  Array<int> &bdr_attr)
622 {
623  Array<int> values_counter;
624  AccumulateAndCountBdrTangentValues(vcoeff, bdr_attr, values_counter);
625 
626  Vector values(Size());
627  for (int i = 0; i < values.Size(); i++)
628  {
629  values(i) = values_counter[i] ? (*this)(i) : 0.0;
630  }
631 
632  // Count the values globally.
633  GroupCommunicator &gcomm = pfes->GroupComm();
634  gcomm.Reduce<int>(values_counter, GroupCommunicator::Sum);
635  // Accumulate the values globally.
636  gcomm.Reduce<double>(values, GroupCommunicator::Sum);
637  // Only the values in the master are guaranteed to be correct!
638  for (int i = 0; i < values.Size(); i++)
639  {
640  if (values_counter[i])
641  {
642  (*this)(i) = values(i)/values_counter[i];
643  }
644  }
645 
646 #ifdef MFEM_DEBUG
647  Array<int> ess_vdofs_marker;
648  pfes->GetEssentialVDofs(bdr_attr, ess_vdofs_marker);
649  for (int i = 0; i < values_counter.Size(); i++)
650  {
651  MFEM_ASSERT(pfes->GetLocalTDofNumber(i) == -1 ||
652  bool(values_counter[i]) == bool(ess_vdofs_marker[i]),
653  "internal error");
654  }
655 #endif
656 }
657 
659  Coefficient *ell_coeff,
660  double Nu,
661  const IntegrationRule *irs[]) const
662 {
663  const_cast<ParGridFunction *>(this)->ExchangeFaceNbrData();
664 
665  int fdof, dim, intorder, k;
666  ElementTransformation *transf;
667  Vector shape, el_dofs, err_val, ell_coeff_val;
668  Array<int> vdofs;
669  IntegrationPoint eip;
670  double error = 0.0;
671 
672  ParMesh *mesh = pfes->GetParMesh();
673  dim = mesh->Dimension();
674 
675  std::map<int,int> local_to_shared;
676  for (int i = 0; i < mesh->GetNSharedFaces(); ++i)
677  {
678  int i_local = mesh->GetSharedFace(i);
679  local_to_shared[i_local] = i;
680  }
681 
682  for (int i = 0; i < mesh->GetNumFaces(); i++)
683  {
684  double shared_face_factor = 1.0;
685  bool shared_face = false;
686  int iel1, iel2, info1, info2;
687  mesh->GetFaceElements(i, &iel1, &iel2);
688  mesh->GetFaceInfos(i, &info1, &info2);
689 
690  intorder = fes->GetFE(iel1)->GetOrder();
691 
692  FaceElementTransformations *face_elem_transf;
693  const FiniteElement *fe1, *fe2;
694  if (info2 >= 0 && iel2 < 0)
695  {
696  int ishared = local_to_shared[i];
697  face_elem_transf = mesh->GetSharedFaceTransformations(ishared);
698  iel2 = face_elem_transf->Elem2No - mesh->GetNE();
699  fe2 = pfes->GetFaceNbrFE(iel2);
700  if ( (k = fe2->GetOrder()) > intorder )
701  {
702  intorder = k;
703  }
704  shared_face = true;
705  shared_face_factor = 0.5;
706  }
707  else
708  {
709  face_elem_transf = mesh->GetFaceElementTransformations(i);
710 
711  if (iel2 >= 0)
712  {
713  fe2 = pfes->GetFE(iel2);
714  if ( (k = fe2->GetOrder()) > intorder )
715  {
716  intorder = k;
717  }
718  }
719  else
720  {
721  fe2 = NULL;
722  }
723  }
724 
725  intorder = 2 * intorder; // <-------------
726  const IntegrationRule *ir;
727  if (irs)
728  {
729  ir = irs[face_elem_transf->GetGeometryType()];
730  }
731  else
732  {
733  ir = &(IntRules.Get(face_elem_transf->GetGeometryType(), intorder));
734  }
735  err_val.SetSize(ir->GetNPoints());
736  ell_coeff_val.SetSize(ir->GetNPoints());
737  // side 1
738  transf = face_elem_transf->Elem1;
739  fe1 = fes->GetFE(iel1);
740  fdof = fe1->GetDof();
741  fes->GetElementVDofs(iel1, vdofs);
742  shape.SetSize(fdof);
743  el_dofs.SetSize(fdof);
744  for (k = 0; k < fdof; k++)
745  if (vdofs[k] >= 0)
746  {
747  el_dofs(k) = (*this)(vdofs[k]);
748  }
749  else
750  {
751  el_dofs(k) = - (*this)(-1-vdofs[k]);
752  }
753  for (int j = 0; j < ir->GetNPoints(); j++)
754  {
755  face_elem_transf->Loc1.Transform(ir->IntPoint(j), eip);
756  fe1->CalcShape(eip, shape);
757  transf->SetIntPoint(&eip);
758  ell_coeff_val(j) = ell_coeff->Eval(*transf, eip);
759  err_val(j) = exsol->Eval(*transf, eip) - (shape * el_dofs);
760  }
761  if (fe2 != NULL)
762  {
763  // side 2
764  transf = face_elem_transf->Elem2;
765  fdof = fe2->GetDof();
766  shape.SetSize(fdof);
767  el_dofs.SetSize(fdof);
768  if (shared_face)
769  {
770  pfes->GetFaceNbrElementVDofs(iel2, vdofs);
771  for (k = 0; k < fdof; k++)
772  if (vdofs[k] >= 0)
773  {
774  el_dofs(k) = face_nbr_data[vdofs[k]];
775  }
776  else
777  {
778  el_dofs(k) = - face_nbr_data[-1-vdofs[k]];
779  }
780  }
781  else
782  {
783  pfes->GetElementVDofs(iel2, vdofs);
784  for (k = 0; k < fdof; k++)
785  if (vdofs[k] >= 0)
786  {
787  el_dofs(k) = (*this)(vdofs[k]);
788  }
789  else
790  {
791  el_dofs(k) = - (*this)(-1 - vdofs[k]);
792  }
793  }
794  for (int j = 0; j < ir->GetNPoints(); j++)
795  {
796  face_elem_transf->Loc2.Transform(ir->IntPoint(j), eip);
797  fe2->CalcShape(eip, shape);
798  transf->SetIntPoint(&eip);
799  ell_coeff_val(j) += ell_coeff->Eval(*transf, eip);
800  ell_coeff_val(j) *= 0.5;
801  err_val(j) -= (exsol->Eval(*transf, eip) - (shape * el_dofs));
802  }
803  }
804  transf = face_elem_transf;
805  for (int j = 0; j < ir->GetNPoints(); j++)
806  {
807  const IntegrationPoint &ip = ir->IntPoint(j);
808  transf->SetIntPoint(&ip);
809  error += shared_face_factor*(ip.weight * Nu * ell_coeff_val(j) *
810  pow(transf->Weight(), 1.0-1.0/(dim-1)) *
811  err_val(j) * err_val(j));
812  }
813  }
814 
815  error = (error < 0.0) ? -sqrt(-error) : sqrt(error);
816  return GlobalLpNorm(2.0, error, pfes->GetComm());
817 }
818 
819 void ParGridFunction::Save(std::ostream &out) const
820 {
821  double *data_ = const_cast<double*>(HostRead());
822  for (int i = 0; i < size; i++)
823  {
824  if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
825  }
826 
827  GridFunction::Save(out);
828 
829  for (int i = 0; i < size; i++)
830  {
831  if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
832  }
833 }
834 
835 #ifdef MFEM_USE_ADIOS2
837  const std::string& variable_name,
838  const adios2stream::data_type type) const
839 {
840  double *data_ = const_cast<double*>(HostRead());
841  for (int i = 0; i < size; i++)
842  {
843  if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
844  }
845 
846  GridFunction::Save(out, variable_name, type);
847 
848  for (int i = 0; i < size; i++)
849  {
850  if (pfes->GetDofSign(i) < 0) { data_[i] = -data_[i]; }
851  }
852 }
853 #endif
854 
855 void ParGridFunction::SaveAsOne(std::ostream &out)
856 {
857  int i, p;
858 
859  MPI_Comm MyComm;
860  MPI_Status status;
861  int MyRank, NRanks;
862 
863  MyComm = pfes -> GetComm();
864 
865  MPI_Comm_size(MyComm, &NRanks);
866  MPI_Comm_rank(MyComm, &MyRank);
867 
868  double **values = new double*[NRanks];
869  int *nv = new int[NRanks];
870  int *nvdofs = new int[NRanks];
871  int *nedofs = new int[NRanks];
872  int *nfdofs = new int[NRanks];
873  int *nrdofs = new int[NRanks];
874 
875  double * h_data = const_cast<double *>(this->HostRead());
876 
877  values[0] = h_data;
878  nv[0] = pfes -> GetVSize();
879  nvdofs[0] = pfes -> GetNVDofs();
880  nedofs[0] = pfes -> GetNEDofs();
881  nfdofs[0] = pfes -> GetNFDofs();
882 
883  if (MyRank == 0)
884  {
885  pfes -> Save(out);
886  out << '\n';
887 
888  for (p = 1; p < NRanks; p++)
889  {
890  MPI_Recv(&nv[p], 1, MPI_INT, p, 455, MyComm, &status);
891  MPI_Recv(&nvdofs[p], 1, MPI_INT, p, 456, MyComm, &status);
892  MPI_Recv(&nedofs[p], 1, MPI_INT, p, 457, MyComm, &status);
893  MPI_Recv(&nfdofs[p], 1, MPI_INT, p, 458, MyComm, &status);
894  values[p] = new double[nv[p]];
895  MPI_Recv(values[p], nv[p], MPI_DOUBLE, p, 460, MyComm, &status);
896  }
897 
898  int vdim = pfes -> GetVDim();
899 
900  for (p = 0; p < NRanks; p++)
901  {
902  nrdofs[p] = nv[p]/vdim - nvdofs[p] - nedofs[p] - nfdofs[p];
903  }
904 
906  {
907  for (int d = 0; d < vdim; d++)
908  {
909  for (p = 0; p < NRanks; p++)
910  for (i = 0; i < nvdofs[p]; i++)
911  {
912  out << *values[p]++ << '\n';
913  }
914 
915  for (p = 0; p < NRanks; p++)
916  for (i = 0; i < nedofs[p]; i++)
917  {
918  out << *values[p]++ << '\n';
919  }
920 
921  for (p = 0; p < NRanks; p++)
922  for (i = 0; i < nfdofs[p]; i++)
923  {
924  out << *values[p]++ << '\n';
925  }
926 
927  for (p = 0; p < NRanks; p++)
928  for (i = 0; i < nrdofs[p]; i++)
929  {
930  out << *values[p]++ << '\n';
931  }
932  }
933  }
934  else
935  {
936  for (p = 0; p < NRanks; p++)
937  for (i = 0; i < nvdofs[p]; i++)
938  for (int d = 0; d < vdim; d++)
939  {
940  out << *values[p]++ << '\n';
941  }
942 
943  for (p = 0; p < NRanks; p++)
944  for (i = 0; i < nedofs[p]; i++)
945  for (int d = 0; d < vdim; d++)
946  {
947  out << *values[p]++ << '\n';
948  }
949 
950  for (p = 0; p < NRanks; p++)
951  for (i = 0; i < nfdofs[p]; i++)
952  for (int d = 0; d < vdim; d++)
953  {
954  out << *values[p]++ << '\n';
955  }
956 
957  for (p = 0; p < NRanks; p++)
958  for (i = 0; i < nrdofs[p]; i++)
959  for (int d = 0; d < vdim; d++)
960  {
961  out << *values[p]++ << '\n';
962  }
963  }
964 
965  for (p = 1; p < NRanks; p++)
966  {
967  values[p] -= nv[p];
968  delete [] values[p];
969  }
970  out.flush();
971  }
972  else
973  {
974  MPI_Send(&nv[0], 1, MPI_INT, 0, 455, MyComm);
975  MPI_Send(&nvdofs[0], 1, MPI_INT, 0, 456, MyComm);
976  MPI_Send(&nedofs[0], 1, MPI_INT, 0, 457, MyComm);
977  MPI_Send(&nfdofs[0], 1, MPI_INT, 0, 458, MyComm);
978  MPI_Send(h_data, nv[0], MPI_DOUBLE, 0, 460, MyComm);
979  }
980 
981  delete [] values;
982  delete [] nv;
983  delete [] nvdofs;
984  delete [] nedofs;
985  delete [] nfdofs;
986  delete [] nrdofs;
987 }
988 
989 double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
990 {
991  double glob_norm;
992 
993  if (p < infinity())
994  {
995  // negative quadrature weights may cause the error to be negative
996  if (loc_norm < 0.0)
997  {
998  loc_norm = -pow(-loc_norm, p);
999  }
1000  else
1001  {
1002  loc_norm = pow(loc_norm, p);
1003  }
1004 
1005  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
1006 
1007  if (glob_norm < 0.0)
1008  {
1009  glob_norm = -pow(-glob_norm, 1.0/p);
1010  }
1011  else
1012  {
1013  glob_norm = pow(glob_norm, 1.0/p);
1014  }
1015  }
1016  else
1017  {
1018  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
1019  }
1020 
1021  return glob_norm;
1022 }
1023 
1025  BilinearFormIntegrator &blfi,
1026  GridFunction &flux, bool wcoef, int subdomain)
1027 {
1028  ParFiniteElementSpace *ffes =
1029  dynamic_cast<ParFiniteElementSpace*>(flux.FESpace());
1030  MFEM_VERIFY(ffes, "the flux FE space must be ParFiniteElementSpace");
1031 
1032  Array<int> count(flux.Size());
1033  SumFluxAndCount(blfi, flux, count, wcoef, subdomain);
1034 
1035  // Accumulate flux and counts in parallel
1036  ffes->GroupComm().Reduce<double>(flux, GroupCommunicator::Sum);
1037  ffes->GroupComm().Bcast<double>(flux);
1038 
1039  ffes->GroupComm().Reduce<int>(count, GroupCommunicator::Sum);
1040  ffes->GroupComm().Bcast<int>(count);
1041 
1042  // complete averaging
1043  for (int i = 0; i < count.Size(); i++)
1044  {
1045  if (count[i] != 0) { flux(i) /= count[i]; }
1046  }
1047 
1048  if (ffes->Nonconforming())
1049  {
1050  // On a partially conforming flux space, project on the conforming space.
1051  // Using this code may lead to worse refinements in ex6, so we do not use
1052  // it by default.
1053 
1054  // Vector conf_flux;
1055  // flux.ConformingProject(conf_flux);
1056  // flux.ConformingProlongate(conf_flux);
1057  }
1058 }
1059 
1060 
1062  const ParGridFunction &x,
1063  ParFiniteElementSpace &smooth_flux_fes,
1064  ParFiniteElementSpace &flux_fes,
1065  Vector &errors,
1066  int norm_p, double solver_tol, int solver_max_it)
1067 {
1068  // Compute fluxes in discontinuous space
1069  GridFunction flux(&flux_fes);
1070  flux = 0.0;
1071 
1072  ParFiniteElementSpace *xfes = x.ParFESpace();
1073  Array<int> xdofs, fdofs;
1074  Vector el_x, el_f;
1075 
1076  for (int i = 0; i < xfes->GetNE(); i++)
1077  {
1078  xfes->GetElementVDofs(i, xdofs);
1079  x.GetSubVector(xdofs, el_x);
1080 
1082  flux_integrator.ComputeElementFlux(*xfes->GetFE(i), *Transf, el_x,
1083  *flux_fes.GetFE(i), el_f, false);
1084 
1085  flux_fes.GetElementVDofs(i, fdofs);
1086  flux.AddElementVector(fdofs, el_f);
1087  }
1088 
1089  // Assemble the linear system for L2 projection into the "smooth" space
1090  ParBilinearForm *a = new ParBilinearForm(&smooth_flux_fes);
1091  ParLinearForm *b = new ParLinearForm(&smooth_flux_fes);
1093 
1094  if (xfes->GetNE())
1095  {
1096  MFEM_VERIFY(smooth_flux_fes.GetFE(0) != NULL,
1097  "Could not obtain FE of smooth flux space.");
1098 
1099  if (smooth_flux_fes.GetFE(0)->GetRangeType() == FiniteElement::SCALAR)
1100  {
1102  vmass->SetVDim(smooth_flux_fes.GetVDim());
1103  a->AddDomainIntegrator(vmass);
1105  }
1106  else
1107  {
1110  }
1111  }
1112 
1113  b->Assemble();
1114  a->Assemble();
1115  a->Finalize();
1116 
1117  // The destination of the projected discontinuous flux
1118  ParGridFunction smooth_flux(&smooth_flux_fes);
1119  smooth_flux = 0.0;
1120 
1121  HypreParMatrix* A = a->ParallelAssemble();
1122  HypreParVector* B = b->ParallelAssemble();
1123  HypreParVector* X = smooth_flux.ParallelProject();
1124 
1125  delete a;
1126  delete b;
1127 
1128  // Define and apply a parallel PCG solver for AX=B with the BoomerAMG
1129  // preconditioner from hypre.
1130  HypreBoomerAMG *amg = new HypreBoomerAMG(*A);
1131  amg->SetPrintLevel(0);
1132  HyprePCG *pcg = new HyprePCG(*A);
1133  pcg->SetTol(solver_tol);
1134  pcg->SetMaxIter(solver_max_it);
1135  pcg->SetPrintLevel(0);
1136  pcg->SetPreconditioner(*amg);
1137  pcg->Mult(*B, *X);
1138 
1139  // Extract the parallel grid function corresponding to the finite element
1140  // approximation X. This is the local solution on each processor.
1141  smooth_flux = *X;
1142 
1143  delete A;
1144  delete B;
1145  delete X;
1146  delete amg;
1147  delete pcg;
1148 
1149  // Proceed through the elements one by one, and find the Lp norm differences
1150  // between the flux as computed per element and the flux projected onto the
1151  // smooth_flux_fes space.
1152  double total_error = 0.0;
1153  errors.SetSize(xfes->GetNE());
1154  for (int i = 0; i < xfes->GetNE(); i++)
1155  {
1156  errors(i) = ComputeElementLpDistance(norm_p, i, smooth_flux, flux);
1157  total_error += pow(errors(i), norm_p);
1158  }
1159 
1160  double glob_error;
1161  MPI_Allreduce(&total_error, &glob_error, 1, MPI_DOUBLE, MPI_SUM,
1162  xfes->GetComm());
1163 
1164  return pow(glob_error, 1.0/norm_p);
1165 }
1166 
1167 } // namespace mfem
1168 
1169 #endif // MFEM_USE_MPI
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for all finite elements.
Definition: fe.hpp:235
int GetNFaceNeighbors() const
Definition: pmesh.hpp:286
void SetTol(double tol)
Definition: hypre.cpp:2619
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:412
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:291
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
Memory< double > data
Definition: vector.hpp:55
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Evaluate the values of all shape functions of a vector finite element in reference space at the given...
Definition: fe.cpp:40
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
Definition: gridfunc.cpp:238
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:198
Geometry::Type GetGeometryType() const
Return the Geometry::Type of the reference element.
Definition: eltrans.hpp:149
Base class for vector Coefficients that optionally depend on time and space.
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
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:2189
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:181
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:414
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Definition: pfespace.cpp:895
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:459
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
double Scale()
Return the scale factor times the optional time dependent function. Returns with when not set by th...
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2387
int GetNSharedFaces() const
Return the number of shared faces (3D), edges (2D), vertices (1D)
Definition: pmesh.cpp:2508
ParFiniteElementSpace * pfes
Points to the same object as fes.
Definition: pgridfunc.hpp:35
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:85
void Copy(Array &copy) const
Create a copy of the internal array to the provided copy.
Definition: array.hpp:831
A specialized ElementTransformation class representing a face and its two neighboring elements...
Definition: eltrans.hpp:467
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:495
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:319
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
Definition: pgridfunc.cpp:855
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
IntegrationPointTransformation Loc2
Definition: eltrans.hpp:511
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:264
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:715
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Definition: pgridfunc.cpp:495
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:87
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, Array< int > &bdr_attr, Array< int > &values_counter)
Definition: gridfunc.cpp:2113
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
Definition: fe.hpp:341
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
Definition: pgridfunc.hpp:44
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
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:2164
ElementTransformation * Elem2
Definition: eltrans.hpp:509
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, double Nu, const IntegrationRule *irs[]=NULL) const
Returns the Face Jumps error for L2 elements.
Definition: pgridfunc.cpp:658
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:969
FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
Definition: mesh.cpp:869
int Size_of_connections() const
Definition: table.hpp:98
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2634
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
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: operator.hpp:92
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
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:427
const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:376
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
Class for parallel linear form.
Definition: plinearform.hpp:26
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:203
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1199
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
ElementTransformation * GetFaceNbrElementTransformation(int i) const
Definition: pfespace.hpp:350
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
IntegrationPointTransformation Loc1
Definition: eltrans.hpp:511
int GetSharedFace(int sface) const
Return the local face index for the given shared face.
Definition: pmesh.cpp:2527
Memory< int > & GetJMemory()
Definition: table.hpp:119
double b
Definition: lissajous.cpp:42
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1119
FaceElementTransformations * GetSharedFaceTransformations(int sf, bool fill2=true)
Definition: pmesh.cpp:2399
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:281
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
void Assemble(int skip_zeros=1)
Assemble the local matrix.
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:123
void AccumulateAndCountZones(Coefficient &coeff, AvgType type, Array< int > &zones_per_vdof)
Accumulates (depending on type) the values of coeff at all shared vdofs and counts in how many zones ...
Definition: gridfunc.cpp:1875
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:188
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true...
Definition: device.hpp:307
int Dimension() const
Definition: mesh.hpp:788
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:587
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2624
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: pgridfunc.cpp:321
int SpaceDimension() const
Definition: mesh.hpp:789
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:460
double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
Definition: pgridfunc.cpp:989
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1031
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1037
PCG solver in hypre.
Definition: hypre.hpp:759
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:582
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 Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
FiniteElementSpace * fes
FE space on which the grid function lives. Owned if fec is not NULL.
Definition: gridfunc.hpp:34
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: pgridfunc.cpp:1024
static void Max(OpData< T >)
Reduce operation Max, instantiated for int and double.
void CalcPhysShape(ElementTransformation &Trans, Vector &shape) const
Evaluate the values of all shape functions of a scalar finite element in physical space at the point ...
Definition: fe.cpp:191
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:298
void Transform(const IntegrationPoint &, IntegrationPoint &)
Definition: eltrans.cpp:533
static bool GetGPUAwareMPI()
Definition: device.hpp:278
FiniteElementCollection * fec
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner()...
Definition: gridfunc.hpp:40
double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, double solver_tol, int solver_max_it)
Definition: pgridfunc.cpp:1061
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:152
double ComputeElementLpDistance(double p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
Definition: gridfunc.cpp:3782
bool Nonconforming() const
Definition: pfespace.hpp:358
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:47
void GetFaceNbrElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: pmesh.cpp:1685
double a
Definition: lissajous.cpp:41
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
Definition: pgridfunc.hpp:39
virtual const char * Name() const
Definition: fe_coll.hpp:61
Class for integration point with weight.
Definition: intrules.hpp:25
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
int dim
Definition: ex24.cpp:53
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
Definition: pgridfunc.cpp:620
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:164
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2639
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Destroy()
Destroy a vector.
Definition: vector.hpp:530
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
ElementTransformation * Elem1
Definition: eltrans.hpp:509
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:138
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:45
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:2252
Class for parallel bilinear form.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:180
int GetFaceNbrVSize() const
Definition: pfespace.hpp:344
double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:384
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:260
Vector data type.
Definition: vector.hpp:51
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:414
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Vector coefficient defined by a vector GridFunction.
int * GetI()
Definition: table.hpp:113
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, bool with_coef=true)
Virtual method required for Zienkiewicz-Zhu type error estimators.
Definition: bilininteg.hpp:214
Class for parallel grid function.
Definition: pgridfunc.hpp:32
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
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2662
int GetRangeType() const
Returns the FiniteElement::RangeType of the element, one of {SCALAR, VECTOR}.
Definition: fe.hpp:332
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:378
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:2293
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:46
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr, Array< int > &values_counter)
Definition: gridfunc.cpp:1971
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:391
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:171
double f(const Vector &p)
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:108