MFEM  v3.4
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, 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 <iostream>
18 #include <limits>
19 using namespace std;
20 
21 namespace mfem
22 {
23 
24 ParGridFunction::ParGridFunction(ParFiniteElementSpace *pf, GridFunction *gf)
25 {
26  fes = pfes = pf;
27  SetDataAndSize(gf->GetData(), gf->Size());
28 }
29 
30 ParGridFunction::ParGridFunction(ParFiniteElementSpace *pf, HypreParVector *tv)
31  : GridFunction(pf), pfes(pf)
32 {
33  Distribute(tv);
34 }
35 
37  const int *partitioning)
38 {
39  const FiniteElementSpace *glob_fes = gf->FESpace();
40  // duplicate the FiniteElementCollection from 'gf'
41  fec = FiniteElementCollection::New(glob_fes->FEColl()->Name());
42  // create a local ParFiniteElementSpace from the global one:
43  fes = pfes = new ParFiniteElementSpace(pmesh, glob_fes, partitioning, fec);
44  SetSize(pfes->GetVSize());
45 
46  if (partitioning)
47  {
48  // Assumption: the map "local element id" -> "global element id" is
49  // increasing, i.e. the local numbering preserves the element order from
50  // the global numbering.
51  Array<int> gvdofs, lvdofs;
52  Vector lnodes;
53  int element_counter = 0;
54  const int MyRank = pfes->GetMyRank();
55  const int glob_ne = glob_fes->GetNE();
56  for (int i = 0; i < glob_ne; i++)
57  {
58  if (partitioning[i] == MyRank)
59  {
60  pfes->GetElementVDofs(element_counter, lvdofs);
61  glob_fes->GetElementVDofs(i, gvdofs);
62  gf->GetSubVector(gvdofs, lnodes);
63  SetSubVector(lvdofs, lnodes);
64  element_counter++;
65  }
66  }
67  }
68 }
69 
70 ParGridFunction::ParGridFunction(ParMesh *pmesh, std::istream &input)
71  : GridFunction(pmesh, input)
72 {
73  // Convert the FiniteElementSpace, fes, to a ParFiniteElementSpace:
74  pfes = new ParFiniteElementSpace(pmesh, fec, fes->GetVDim(),
75  fes->GetOrdering());
76  delete fes;
77  fes = pfes;
78 }
79 
81 {
84 }
85 
87 {
90  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
91  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
92 }
93 
95 {
98  pfes = f;
99 }
100 
102 {
104  GridFunction::MakeRef(f, v);
105  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
106  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
107 }
108 
110 {
112  GridFunction::MakeRef(f, v);
113  pfes = f;
114 }
115 
117 {
119  GridFunction::MakeRef(f, v, v_offset);
120  pfes = dynamic_cast<ParFiniteElementSpace*>(f);
121  MFEM_ASSERT(pfes != NULL, "not a ParFiniteElementSpace");
122 }
123 
125 {
127  GridFunction::MakeRef(f, v, v_offset);
128  pfes = f;
129 }
130 
132 {
133  pfes->GetProlongationMatrix()->Mult(*tv, *this);
134 }
135 
136 void ParGridFunction::AddDistribute(double a, const Vector *tv)
137 {
138  pfes->Dof_TrueDof_Matrix()->Mult(a, *tv, 1.0, *this);
139 }
140 
142 {
144  GetTrueDofs(*tv);
145  return tv;
146 }
147 
149 {
151  pfes->DivideByGroupSize(tv);
152 }
153 
155 {
157  pfes->DivideByGroupSize(tv);
158 }
159 
161 {
163  ParallelAverage(*tv);
164  return tv;
165 }
166 
168 {
169  pfes->GetRestrictionMatrix()->Mult(*this, tv);
170 }
171 
173 {
174  pfes->GetRestrictionMatrix()->Mult(*this, tv);
175 }
176 
178 {
180  ParallelProject(*tv);
181  return tv;
182 }
183 
185 {
187 }
188 
190 {
192 }
193 
195 {
197  ParallelAssemble(*tv);
198  return tv;
199 }
200 
202 {
204 
205  if (pfes->GetFaceNbrVSize() <= 0)
206  {
207  return;
208  }
209 
210  ParMesh *pmesh = pfes->GetParMesh();
211 
214 
215  int *send_offset = pfes->send_face_nbr_ldof.GetI();
216  int *send_ldof = pfes->send_face_nbr_ldof.GetJ();
217  int *recv_offset = pfes->face_nbr_ldof.GetI();
218  MPI_Comm MyComm = pfes->GetComm();
219 
220  int num_face_nbrs = pmesh->GetNFaceNeighbors();
221  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
222  MPI_Request *send_requests = requests;
223  MPI_Request *recv_requests = requests + num_face_nbrs;
224  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
225 
226  for (int i = 0; i < send_data.Size(); i++)
227  {
228  send_data[i] = data[send_ldof[i]];
229  }
230 
231  for (int fn = 0; fn < num_face_nbrs; fn++)
232  {
233  int nbr_rank = pmesh->GetFaceNbrRank(fn);
234  int tag = 0;
235 
236  MPI_Isend(&send_data(send_offset[fn]),
237  send_offset[fn+1] - send_offset[fn],
238  MPI_DOUBLE, nbr_rank, tag, MyComm, &send_requests[fn]);
239 
240  MPI_Irecv(&face_nbr_data(recv_offset[fn]),
241  recv_offset[fn+1] - recv_offset[fn],
242  MPI_DOUBLE, nbr_rank, tag, MyComm, &recv_requests[fn]);
243  }
244 
245  MPI_Waitall(num_face_nbrs, send_requests, statuses);
246  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
247 
248  delete [] statuses;
249  delete [] requests;
250 }
251 
252 double ParGridFunction::GetValue(int i, const IntegrationPoint &ip, int vdim)
253 const
254 {
255  Array<int> dofs;
256  Vector DofVal, LocVec;
257  int nbr_el_no = i - pfes->GetParMesh()->GetNE();
258  if (nbr_el_no >= 0)
259  {
260  int fes_vdim = pfes->GetVDim();
261  pfes->GetFaceNbrElementVDofs(nbr_el_no, dofs);
262  if (fes_vdim > 1)
263  {
264  int s = dofs.Size()/fes_vdim;
265  Array<int> _dofs(&dofs[(vdim-1)*s], s);
266  face_nbr_data.GetSubVector(_dofs, LocVec);
267  DofVal.SetSize(s);
268  }
269  else
270  {
271  face_nbr_data.GetSubVector(dofs, LocVec);
272  DofVal.SetSize(dofs.Size());
273  }
274  pfes->GetFaceNbrFE(nbr_el_no)->CalcShape(ip, DofVal);
275  }
276  else
277  {
278  fes->GetElementDofs(i, dofs);
279  fes->DofsToVDofs(vdim-1, dofs);
280  DofVal.SetSize(dofs.Size());
281  const FiniteElement *fe = fes->GetFE(i);
282  MFEM_ASSERT(fe->GetMapType() == FiniteElement::VALUE, "invalid FE map type");
283  fe->CalcShape(ip, DofVal);
284  GetSubVector(dofs, LocVec);
285  }
286 
287  return (DofVal * LocVec);
288 }
289 
291 {
292  DeltaCoefficient *delta_c = dynamic_cast<DeltaCoefficient *>(&coeff);
293 
294  if (delta_c == NULL)
295  {
297  }
298  else
299  {
300  double loc_integral, glob_integral;
301 
302  ProjectDeltaCoefficient(*delta_c, loc_integral);
303 
304  MPI_Allreduce(&loc_integral, &glob_integral, 1, MPI_DOUBLE, MPI_SUM,
305  pfes->GetComm());
306 
307  (*this) *= (delta_c->Scale() / glob_integral);
308  }
309 }
310 
312 {
313  // local maximal element attribute for each dof
314  Array<int> ldof_attr;
315 
316  // local projection
317  GridFunction::ProjectDiscCoefficient(coeff, ldof_attr);
318 
319  // global maximal element attribute for each dof
320  Array<int> gdof_attr;
321  ldof_attr.Copy(gdof_attr);
322  GroupCommunicator &gcomm = pfes->GroupComm();
323  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Max);
324  gcomm.Bcast(gdof_attr);
325 
326  // set local value to zero if global maximal element attribute is larger than
327  // the local one, and mark (in gdof_attr) if we have the correct value
328  for (int i = 0; i < pfes->GetVSize(); i++)
329  {
330  if (gdof_attr[i] > ldof_attr[i])
331  {
332  (*this)(i) = 0.0;
333  gdof_attr[i] = 0;
334  }
335  else
336  {
337  gdof_attr[i] = 1;
338  }
339  }
340 
341  // parallel averaging plus interpolation to determine final values
343  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Sum);
344  gcomm.Bcast(gdof_attr);
345  for (int i = 0; i < fes->GetVSize(); i++)
346  {
347  (*this)(i) /= gdof_attr[i];
348  }
349  this->ParallelAssemble(*tv);
350  this->Distribute(tv);
351  delete tv;
352 }
353 
354 
356 {
357  // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
358  // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
359 
360  // Number of zones that contain a given dof.
361  Array<int> zones_per_vdof;
362  AccumulateAndCountZones(coeff, type, zones_per_vdof);
363 
364  // Count the zones globally.
365  GroupCommunicator &gcomm = pfes->GroupComm();
366  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
367  gcomm.Bcast(zones_per_vdof);
368  // Accumulate for all tdofs.
369  HypreParVector *tv = this->ParallelAssemble();
370  this->Distribute(tv);
371  delete tv;
372 
373  ComputeMeans(type, zones_per_vdof);
374 }
375 
377  AvgType type)
378 {
379  // Harmonic (x1 ... xn) = [ (1/x1 + ... + 1/xn) / n ]^-1.
380  // Arithmetic(x1 ... xn) = (x1 + ... + xn) / n.
381 
382  // Number of zones that contain a given dof.
383  Array<int> zones_per_vdof;
384  AccumulateAndCountZones(vcoeff, type, zones_per_vdof);
385 
386  // Count the zones globally.
387  GroupCommunicator &gcomm = pfes->GroupComm();
388  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
389  gcomm.Bcast(zones_per_vdof);
390  // Accumulate for all tdofs.
391  HypreParVector *tv = this->ParallelAssemble();
392  this->Distribute(tv);
393  delete tv;
394 
395  ComputeMeans(type, zones_per_vdof);
396 }
397 
398 void ParGridFunction::Save(std::ostream &out) const
399 {
400  for (int i = 0; i < size; i++)
401  {
402  if (pfes->GetDofSign(i) < 0) { data[i] = -data[i]; }
403  }
404 
405  GridFunction::Save(out);
406 
407  for (int i = 0; i < size; i++)
408  {
409  if (pfes->GetDofSign(i) < 0) { data[i] = -data[i]; }
410  }
411 }
412 
413 void ParGridFunction::SaveAsOne(std::ostream &out)
414 {
415  int i, p;
416 
417  MPI_Comm MyComm;
418  MPI_Status status;
419  int MyRank, NRanks;
420 
421  MyComm = pfes -> GetComm();
422 
423  MPI_Comm_size(MyComm, &NRanks);
424  MPI_Comm_rank(MyComm, &MyRank);
425 
426  double **values = new double*[NRanks];
427  int *nv = new int[NRanks];
428  int *nvdofs = new int[NRanks];
429  int *nedofs = new int[NRanks];
430  int *nfdofs = new int[NRanks];
431  int *nrdofs = new int[NRanks];
432 
433  values[0] = data;
434  nv[0] = pfes -> GetVSize();
435  nvdofs[0] = pfes -> GetNVDofs();
436  nedofs[0] = pfes -> GetNEDofs();
437  nfdofs[0] = pfes -> GetNFDofs();
438 
439  if (MyRank == 0)
440  {
441  pfes -> Save(out);
442  out << '\n';
443 
444  for (p = 1; p < NRanks; p++)
445  {
446  MPI_Recv(&nv[p], 1, MPI_INT, p, 455, MyComm, &status);
447  MPI_Recv(&nvdofs[p], 1, MPI_INT, p, 456, MyComm, &status);
448  MPI_Recv(&nedofs[p], 1, MPI_INT, p, 457, MyComm, &status);
449  MPI_Recv(&nfdofs[p], 1, MPI_INT, p, 458, MyComm, &status);
450  values[p] = new double[nv[p]];
451  MPI_Recv(values[p], nv[p], MPI_DOUBLE, p, 460, MyComm, &status);
452  }
453 
454  int vdim = pfes -> GetVDim();
455 
456  for (p = 0; p < NRanks; p++)
457  {
458  nrdofs[p] = nv[p]/vdim - nvdofs[p] - nedofs[p] - nfdofs[p];
459  }
460 
462  {
463  for (int d = 0; d < vdim; d++)
464  {
465  for (p = 0; p < NRanks; p++)
466  for (i = 0; i < nvdofs[p]; i++)
467  {
468  out << *values[p]++ << '\n';
469  }
470 
471  for (p = 0; p < NRanks; p++)
472  for (i = 0; i < nedofs[p]; i++)
473  {
474  out << *values[p]++ << '\n';
475  }
476 
477  for (p = 0; p < NRanks; p++)
478  for (i = 0; i < nfdofs[p]; i++)
479  {
480  out << *values[p]++ << '\n';
481  }
482 
483  for (p = 0; p < NRanks; p++)
484  for (i = 0; i < nrdofs[p]; i++)
485  {
486  out << *values[p]++ << '\n';
487  }
488  }
489  }
490  else
491  {
492  for (p = 0; p < NRanks; p++)
493  for (i = 0; i < nvdofs[p]; i++)
494  for (int d = 0; d < vdim; d++)
495  {
496  out << *values[p]++ << '\n';
497  }
498 
499  for (p = 0; p < NRanks; p++)
500  for (i = 0; i < nedofs[p]; i++)
501  for (int d = 0; d < vdim; d++)
502  {
503  out << *values[p]++ << '\n';
504  }
505 
506  for (p = 0; p < NRanks; p++)
507  for (i = 0; i < nfdofs[p]; i++)
508  for (int d = 0; d < vdim; d++)
509  {
510  out << *values[p]++ << '\n';
511  }
512 
513  for (p = 0; p < NRanks; p++)
514  for (i = 0; i < nrdofs[p]; i++)
515  for (int d = 0; d < vdim; d++)
516  {
517  out << *values[p]++ << '\n';
518  }
519  }
520 
521  for (p = 1; p < NRanks; p++)
522  {
523  values[p] -= nv[p];
524  delete [] values[p];
525  }
526  out.flush();
527  }
528  else
529  {
530  MPI_Send(&nv[0], 1, MPI_INT, 0, 455, MyComm);
531  MPI_Send(&nvdofs[0], 1, MPI_INT, 0, 456, MyComm);
532  MPI_Send(&nedofs[0], 1, MPI_INT, 0, 457, MyComm);
533  MPI_Send(&nfdofs[0], 1, MPI_INT, 0, 458, MyComm);
534  MPI_Send(data, nv[0], MPI_DOUBLE, 0, 460, MyComm);
535  }
536 
537  delete [] values;
538  delete [] nv;
539  delete [] nvdofs;
540  delete [] nedofs;
541  delete [] nfdofs;
542  delete [] nrdofs;
543 }
544 
545 double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
546 {
547  double glob_norm;
548 
549  if (p < infinity())
550  {
551  // negative quadrature weights may cause the error to be negative
552  if (loc_norm < 0.0)
553  {
554  loc_norm = -pow(-loc_norm, p);
555  }
556  else
557  {
558  loc_norm = pow(loc_norm, p);
559  }
560 
561  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
562 
563  if (glob_norm < 0.0)
564  {
565  glob_norm = -pow(-glob_norm, 1.0/p);
566  }
567  else
568  {
569  glob_norm = pow(glob_norm, 1.0/p);
570  }
571  }
572  else
573  {
574  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
575  }
576 
577  return glob_norm;
578 }
579 
580 
583  GridFunction &flux, int wcoef, int subdomain)
584 {
585  ParFiniteElementSpace *ffes =
586  dynamic_cast<ParFiniteElementSpace*>(flux.FESpace());
587  MFEM_VERIFY(ffes, "the flux FE space must be ParFiniteElementSpace");
588 
589  Array<int> count(flux.Size());
590  SumFluxAndCount(blfi, flux, count, wcoef, subdomain);
591 
592  if (ffes->Conforming()) // FIXME: nonconforming
593  {
594  // Accumulate flux and counts in parallel
595 
596  ffes->GroupComm().Reduce<double>(flux, GroupCommunicator::Sum);
597  ffes->GroupComm().Bcast<double>(flux);
598 
599  ffes->GroupComm().Reduce<int>(count, GroupCommunicator::Sum);
600  ffes->GroupComm().Bcast<int>(count);
601  }
602  else
603  {
604  MFEM_ABORT("Averaging on processor boundaries not implemented for "
605  "NC meshes yet.\n"
606  "Use L2ZZErrorEstimator() instead of ZZErrorEstimator().");
607  }
608 
609  // complete averaging
610  for (int i = 0; i < count.Size(); i++)
611  {
612  if (count[i] != 0) { flux(i) /= count[i]; }
613  }
614 
615  if (ffes->Nonconforming())
616  {
617  // On a partially conforming flux space, project on the conforming space.
618  // Using this code may lead to worse refinements in ex6, so we do not use
619  // it by default.
620 
621  // Vector conf_flux;
622  // flux.ConformingProject(conf_flux);
623  // flux.ConformingProlongate(conf_flux);
624  }
625 }
626 
627 
629  const ParGridFunction &x,
630  ParFiniteElementSpace &smooth_flux_fes,
631  ParFiniteElementSpace &flux_fes,
632  Vector &errors,
633  int norm_p, double solver_tol, int solver_max_it)
634 {
635  // Compute fluxes in discontinuous space
636  GridFunction flux(&flux_fes);
637  flux = 0.0;
638 
639  ParFiniteElementSpace *xfes = x.ParFESpace();
640  Array<int> xdofs, fdofs;
641  Vector el_x, el_f;
642 
643  for (int i = 0; i < xfes->GetNE(); i++)
644  {
645  xfes->GetElementVDofs(i, xdofs);
646  x.GetSubVector(xdofs, el_x);
647 
649  flux_integrator.ComputeElementFlux(*xfes->GetFE(i), *Transf, el_x,
650  *flux_fes.GetFE(i), el_f, false);
651 
652  flux_fes.GetElementVDofs(i, fdofs);
653  flux.AddElementVector(fdofs, el_f);
654  }
655 
656  // Assemble the linear system for L2 projection into the "smooth" space
657  ParBilinearForm *a = new ParBilinearForm(&smooth_flux_fes);
658  ParLinearForm *b = new ParLinearForm(&smooth_flux_fes);
660 
661  if (smooth_flux_fes.GetFE(0)->GetRangeType() == FiniteElement::SCALAR)
662  {
665  }
666  else
667  {
670  }
671 
672  b->Assemble();
673  a->Assemble();
674  a->Finalize();
675 
676  // The destination of the projected discontinuous flux
677  ParGridFunction smooth_flux(&smooth_flux_fes);
678  smooth_flux = 0.0;
679 
682  HypreParVector* X = smooth_flux.ParallelProject();
683 
684  delete a;
685  delete b;
686 
687  // Define and apply a parallel PCG solver for AX=B with the BoomerAMG
688  // preconditioner from hypre.
689  HypreBoomerAMG *amg = new HypreBoomerAMG(*A);
690  amg->SetPrintLevel(0);
691  HyprePCG *pcg = new HyprePCG(*A);
692  pcg->SetTol(solver_tol);
693  pcg->SetMaxIter(solver_max_it);
694  pcg->SetPrintLevel(0);
695  pcg->SetPreconditioner(*amg);
696  pcg->Mult(*B, *X);
697 
698  // Extract the parallel grid function corresponding to the finite element
699  // approximation X. This is the local solution on each processor.
700  smooth_flux = *X;
701 
702  delete A;
703  delete B;
704  delete X;
705  delete amg;
706  delete pcg;
707 
708  // Proceed through the elements one by one, and find the Lp norm differences
709  // between the flux as computed per element and the flux projected onto the
710  // smooth_flux_fes space.
711  double total_error = 0.0;
712  errors.SetSize(xfes->GetNE());
713  for (int i = 0; i < xfes->GetNE(); i++)
714  {
715  errors(i) = ComputeElementLpDistance(norm_p, i, smooth_flux, flux);
716  total_error += pow(errors(i), norm_p);
717  }
718 
719  double glob_error;
720  MPI_Allreduce(&total_error, &glob_error, 1, MPI_DOUBLE, MPI_SUM,
721  xfes->GetComm());
722 
723  return pow(glob_error, 1.0/norm_p);
724 }
725 
726 }
727 
728 #endif // MFEM_USE_MPI
Abstract class for Finite Elements.
Definition: fe.hpp:140
int GetNFaceNeighbors() const
Definition: pmesh.hpp:161
void SetTol(double tol)
Definition: hypre.cpp:2170
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:265
int Size() const
Logical size of the array.
Definition: array.hpp:133
For scalar fields; preserves point values.
Definition: fe.hpp:180
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:494
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:253
int * GetJ()
Definition: table.hpp:114
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:277
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:194
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1314
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:177
virtual const Operator * GetProlongationMatrix() const
Definition: pfespace.cpp:742
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
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:171
double Scale()
Return the scale set by SetScale() multiplied by the time-dependent function specified by SetFunction...
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1508
ParFiniteElementSpace * pfes
Definition: pgridfunc.hpp:35
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:190
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:458
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
Definition: pgridfunc.cpp:413
Delta function coefficient.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:252
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:398
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:557
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:311
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:86
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int GetMapType() const
Definition: fe.hpp:237
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:290
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:1291
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:129
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:1005
int Size_of_connections() const
Definition: table.hpp:98
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2185
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:52
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:974
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:101
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2349
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:277
The BoomerAMG solver in hypre.
Definition: hypre.hpp:796
Class for parallel linear form.
Definition: plinearform.hpp:26
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:229
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:1008
MPI_Comm GetComm() const
Definition: pfespace.hpp:235
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:80
void SetPrintLevel(int print_level)
Definition: hypre.hpp:837
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:267
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:1233
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void ComputeElementFlux(const FiniteElement &el, ElementTransformation &Trans, Vector &u, const FiniteElement &fluxelem, Vector &flux, int with_coef=1)
Definition: bilininteg.hpp:72
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:1201
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:182
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:242
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:337
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:546
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2175
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:301
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:545
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
PCG solver in hypre.
Definition: hypre.hpp:656
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:486
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
FiniteElementSpace * fes
FE space on which grid function lives.
Definition: gridfunc.hpp:31
static void Max(OpData< T >)
Reduce operation Max, instantiated for int and double.
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:284
FiniteElementCollection * fec
Used when the grid function is read from a file.
Definition: gridfunc.hpp:34
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:628
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:131
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:147
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:2707
bool Nonconforming() const
Definition: pfespace.hpp:339
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:44
virtual const char * Name() const
Definition: fe_coll.hpp:50
Class for integration point with weight.
Definition: intrules.hpp:25
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:322
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:160
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2190
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void Destroy()
Destroy a vector.
Definition: vector.hpp:347
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1335
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:141
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:136
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1377
Class for parallel bilinear form.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:174
int GetFaceNbrVSize() const
Definition: pfespace.hpp:327
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:231
Vector data type.
Definition: vector.hpp:48
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:267
Vector coefficient defined by a vector GridFunction.
int * GetI()
Definition: table.hpp:113
Class for parallel grid function.
Definition: pgridfunc.hpp:32
double * data
Definition: vector.hpp:53
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:64
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2211
int GetRangeType() const
Definition: fe.hpp:233
Class for parallel meshes.
Definition: pmesh.hpp:32
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:1839
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:34
Integrator for (Q u, v) for VectorFiniteElements.
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: pgridfunc.cpp:581
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:167
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:85
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:106