MFEM  v3.3
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  int * partitioning)
38 {
39  // duplicate the FiniteElementCollection from 'gf'
41  fes = pfes = new ParFiniteElementSpace(pmesh, fec, gf->FESpace()->GetVDim(),
42  gf->FESpace()->GetOrdering());
43  SetSize(pfes->GetVSize());
44 
45  if (partitioning)
46  {
47  Array<int> gvdofs, lvdofs;
48  Vector lnodes;
49  int element_counter = 0;
50  Mesh & mesh(*gf->FESpace()->GetMesh());
51  int MyRank;
52  MPI_Comm_rank(pfes->GetComm(), &MyRank);
53  for (int i = 0; i < mesh.GetNE(); i++)
54  if (partitioning[i] == MyRank)
55  {
56  pfes->GetElementVDofs(element_counter, lvdofs);
57  gf->FESpace()->GetElementVDofs(i, gvdofs);
58  gf->GetSubVector(gvdofs, lnodes);
59  SetSubVector(lvdofs, lnodes);
60  element_counter++;
61  }
62  }
63 }
64 
66 {
69 }
70 
72 {
75  pfes = f;
76 }
77 
79 {
82  pfes = f;
83 }
84 
86 {
88  GridFunction::MakeRef(f, v, v_offset);
89  pfes = f;
90 }
91 
93 {
94  pfes->Dof_TrueDof_Matrix()->Mult(*tv, *this);
95 }
96 
97 void ParGridFunction::AddDistribute(double a, const Vector *tv)
98 {
99  pfes->Dof_TrueDof_Matrix()->Mult(a, *tv, 1.0, *this);
100 }
101 
103 {
105  GetTrueDofs(*tv);
106  return tv;
107 }
108 
110 {
111  pfes->Dof_TrueDof_Matrix()->MultTranspose(*this, tv);
112  pfes->DivideByGroupSize(tv);
113 }
114 
116 {
117  pfes->Dof_TrueDof_Matrix()->MultTranspose(*this, tv);
118  pfes->DivideByGroupSize(tv);
119 }
120 
122 {
124  ParallelAverage(*tv);
125  return tv;
126 }
127 
129 {
130  pfes->GetRestrictionMatrix()->Mult(*this, tv);
131 }
132 
134 {
135  pfes->GetRestrictionMatrix()->Mult(*this, tv);
136 }
137 
139 {
141  ParallelProject(*tv);
142  return tv;
143 }
144 
146 {
147  pfes->Dof_TrueDof_Matrix()->MultTranspose(*this, tv);
148 }
149 
151 {
152  pfes->Dof_TrueDof_Matrix()->MultTranspose(*this, tv);
153 }
154 
156 {
158  ParallelAssemble(*tv);
159  return tv;
160 }
161 
163 {
165 
166  if (pfes->GetFaceNbrVSize() <= 0)
167  {
168  return;
169  }
170 
171  ParMesh *pmesh = pfes->GetParMesh();
172 
175 
176  int *send_offset = pfes->send_face_nbr_ldof.GetI();
177  int *send_ldof = pfes->send_face_nbr_ldof.GetJ();
178  int *recv_offset = pfes->face_nbr_ldof.GetI();
179  MPI_Comm MyComm = pfes->GetComm();
180 
181  int num_face_nbrs = pmesh->GetNFaceNeighbors();
182  MPI_Request *requests = new MPI_Request[2*num_face_nbrs];
183  MPI_Request *send_requests = requests;
184  MPI_Request *recv_requests = requests + num_face_nbrs;
185  MPI_Status *statuses = new MPI_Status[num_face_nbrs];
186 
187  for (int i = 0; i < send_data.Size(); i++)
188  {
189  send_data[i] = data[send_ldof[i]];
190  }
191 
192  for (int fn = 0; fn < num_face_nbrs; fn++)
193  {
194  int nbr_rank = pmesh->GetFaceNbrRank(fn);
195  int tag = 0;
196 
197  MPI_Isend(&send_data(send_offset[fn]),
198  send_offset[fn+1] - send_offset[fn],
199  MPI_DOUBLE, nbr_rank, tag, MyComm, &send_requests[fn]);
200 
201  MPI_Irecv(&face_nbr_data(recv_offset[fn]),
202  recv_offset[fn+1] - recv_offset[fn],
203  MPI_DOUBLE, nbr_rank, tag, MyComm, &recv_requests[fn]);
204  }
205 
206  MPI_Waitall(num_face_nbrs, send_requests, statuses);
207  MPI_Waitall(num_face_nbrs, recv_requests, statuses);
208 
209  delete [] statuses;
210  delete [] requests;
211 }
212 
213 double ParGridFunction::GetValue(int i, const IntegrationPoint &ip, int vdim)
214 const
215 {
216  Array<int> dofs;
217  Vector DofVal, LocVec;
218  int nbr_el_no = i - pfes->GetParMesh()->GetNE();
219  if (nbr_el_no >= 0)
220  {
221  int fes_vdim = pfes->GetVDim();
222  pfes->GetFaceNbrElementVDofs(nbr_el_no, dofs);
223  if (fes_vdim > 1)
224  {
225  int s = dofs.Size()/fes_vdim;
226  Array<int> _dofs(&dofs[(vdim-1)*s], s);
227  face_nbr_data.GetSubVector(_dofs, LocVec);
228  DofVal.SetSize(s);
229  }
230  else
231  {
232  face_nbr_data.GetSubVector(dofs, LocVec);
233  DofVal.SetSize(dofs.Size());
234  }
235  pfes->GetFaceNbrFE(nbr_el_no)->CalcShape(ip, DofVal);
236  }
237  else
238  {
239  fes->GetElementDofs(i, dofs);
240  fes->DofsToVDofs(vdim-1, dofs);
241  DofVal.SetSize(dofs.Size());
242  fes->GetFE(i)->CalcShape(ip, DofVal);
243  GetSubVector(dofs, LocVec);
244  }
245 
246  return (DofVal * LocVec);
247 }
248 
250 {
251  DeltaCoefficient *delta_c = dynamic_cast<DeltaCoefficient *>(&coeff);
252 
253  if (delta_c == NULL)
254  {
256  }
257  else
258  {
259  double loc_integral, glob_integral;
260 
261  ProjectDeltaCoefficient(*delta_c, loc_integral);
262 
263  MPI_Allreduce(&loc_integral, &glob_integral, 1, MPI_DOUBLE, MPI_SUM,
264  pfes->GetComm());
265 
266  (*this) *= (delta_c->Scale() / glob_integral);
267  }
268 }
269 
271 {
272  // local maximal element attribute for each dof
273  Array<int> ldof_attr;
274 
275  // local projection
276  GridFunction::ProjectDiscCoefficient(coeff, ldof_attr);
277 
278  // global maximal element attribute for each dof
279  Array<int> gdof_attr;
280  ldof_attr.Copy(gdof_attr);
281  GroupCommunicator &gcomm = pfes->GroupComm();
282  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Max);
283  gcomm.Bcast(gdof_attr);
284 
285  // set local value to zero if global maximal element attribute is larger than
286  // the local one, and mark (in gdof_attr) if we have the correct value
287  for (int i = 0; i < pfes->GetVSize(); i++)
288  {
289  if (gdof_attr[i] > ldof_attr[i])
290  {
291  (*this)(i) = 0.0;
292  gdof_attr[i] = 0;
293  }
294  else
295  {
296  gdof_attr[i] = 1;
297  }
298  }
299 
300  // parallel averaging plus interpolation to determine final values
302  gcomm.Reduce<int>(gdof_attr, GroupCommunicator::Sum);
303  gcomm.Bcast(gdof_attr);
304  for (int i = 0; i < fes->GetVSize(); i++)
305  {
306  (*this)(i) /= gdof_attr[i];
307  }
308  this->ParallelAssemble(*tv);
309  this->Distribute(tv);
310  delete tv;
311 }
312 
313 
314 void ParGridFunction::Save(std::ostream &out) const
315 {
316  for (int i = 0; i < size; i++)
317  {
318  if (pfes->GetDofSign(i) < 0) { data[i] = -data[i]; }
319  }
320 
321  GridFunction::Save(out);
322 
323  for (int i = 0; i < size; i++)
324  {
325  if (pfes->GetDofSign(i) < 0) { data[i] = -data[i]; }
326  }
327 }
328 
329 void ParGridFunction::SaveAsOne(std::ostream &out)
330 {
331  int i, p;
332 
333  MPI_Comm MyComm;
334  MPI_Status status;
335  int MyRank, NRanks;
336 
337  MyComm = pfes -> GetComm();
338 
339  MPI_Comm_size(MyComm, &NRanks);
340  MPI_Comm_rank(MyComm, &MyRank);
341 
342  double **values = new double*[NRanks];
343  int *nv = new int[NRanks];
344  int *nvdofs = new int[NRanks];
345  int *nedofs = new int[NRanks];
346  int *nfdofs = new int[NRanks];
347  int *nrdofs = new int[NRanks];
348 
349  values[0] = data;
350  nv[0] = pfes -> GetVSize();
351  nvdofs[0] = pfes -> GetNVDofs();
352  nedofs[0] = pfes -> GetNEDofs();
353  nfdofs[0] = pfes -> GetNFDofs();
354 
355  if (MyRank == 0)
356  {
357  pfes -> Save(out);
358  out << '\n';
359 
360  for (p = 1; p < NRanks; p++)
361  {
362  MPI_Recv(&nv[p], 1, MPI_INT, p, 455, MyComm, &status);
363  MPI_Recv(&nvdofs[p], 1, MPI_INT, p, 456, MyComm, &status);
364  MPI_Recv(&nedofs[p], 1, MPI_INT, p, 457, MyComm, &status);
365  MPI_Recv(&nfdofs[p], 1, MPI_INT, p, 458, MyComm, &status);
366  values[p] = new double[nv[p]];
367  MPI_Recv(values[p], nv[p], MPI_DOUBLE, p, 460, MyComm, &status);
368  }
369 
370  int vdim = pfes -> GetVDim();
371 
372  for (p = 0; p < NRanks; p++)
373  {
374  nrdofs[p] = nv[p]/vdim - nvdofs[p] - nedofs[p] - nfdofs[p];
375  }
376 
378  {
379  for (int d = 0; d < vdim; d++)
380  {
381  for (p = 0; p < NRanks; p++)
382  for (i = 0; i < nvdofs[p]; i++)
383  {
384  out << *values[p]++ << '\n';
385  }
386 
387  for (p = 0; p < NRanks; p++)
388  for (i = 0; i < nedofs[p]; i++)
389  {
390  out << *values[p]++ << '\n';
391  }
392 
393  for (p = 0; p < NRanks; p++)
394  for (i = 0; i < nfdofs[p]; i++)
395  {
396  out << *values[p]++ << '\n';
397  }
398 
399  for (p = 0; p < NRanks; p++)
400  for (i = 0; i < nrdofs[p]; i++)
401  {
402  out << *values[p]++ << '\n';
403  }
404  }
405  }
406  else
407  {
408  for (p = 0; p < NRanks; p++)
409  for (i = 0; i < nvdofs[p]; i++)
410  for (int d = 0; d < vdim; d++)
411  {
412  out << *values[p]++ << '\n';
413  }
414 
415  for (p = 0; p < NRanks; p++)
416  for (i = 0; i < nedofs[p]; i++)
417  for (int d = 0; d < vdim; d++)
418  {
419  out << *values[p]++ << '\n';
420  }
421 
422  for (p = 0; p < NRanks; p++)
423  for (i = 0; i < nfdofs[p]; i++)
424  for (int d = 0; d < vdim; d++)
425  {
426  out << *values[p]++ << '\n';
427  }
428 
429  for (p = 0; p < NRanks; p++)
430  for (i = 0; i < nrdofs[p]; i++)
431  for (int d = 0; d < vdim; d++)
432  {
433  out << *values[p]++ << '\n';
434  }
435  }
436 
437  for (p = 1; p < NRanks; p++)
438  {
439  values[p] -= nv[p];
440  delete [] values[p];
441  }
442  out.flush();
443  }
444  else
445  {
446  MPI_Send(&nv[0], 1, MPI_INT, 0, 455, MyComm);
447  MPI_Send(&nvdofs[0], 1, MPI_INT, 0, 456, MyComm);
448  MPI_Send(&nedofs[0], 1, MPI_INT, 0, 457, MyComm);
449  MPI_Send(&nfdofs[0], 1, MPI_INT, 0, 458, MyComm);
450  MPI_Send(data, nv[0], MPI_DOUBLE, 0, 460, MyComm);
451  }
452 
453  delete [] values;
454  delete [] nv;
455  delete [] nvdofs;
456  delete [] nedofs;
457  delete [] nfdofs;
458  delete [] nrdofs;
459 }
460 
461 double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
462 {
463  double glob_norm;
464 
465  if (p < numeric_limits<double>::infinity())
466  {
467  // negative quadrature weights may cause the error to be negative
468  if (loc_norm < 0.0)
469  {
470  loc_norm = -pow(-loc_norm, p);
471  }
472  else
473  {
474  loc_norm = pow(loc_norm, p);
475  }
476 
477  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_SUM, comm);
478 
479  if (glob_norm < 0.0)
480  {
481  glob_norm = -pow(-glob_norm, 1.0/p);
482  }
483  else
484  {
485  glob_norm = pow(glob_norm, 1.0/p);
486  }
487  }
488  else
489  {
490  MPI_Allreduce(&loc_norm, &glob_norm, 1, MPI_DOUBLE, MPI_MAX, comm);
491  }
492 
493  return glob_norm;
494 }
495 
496 
499  GridFunction &flux, int wcoef, int subdomain)
500 {
501  ParFiniteElementSpace *ffes =
502  dynamic_cast<ParFiniteElementSpace*>(flux.FESpace());
503  MFEM_VERIFY(ffes, "the flux FE space must be ParFiniteElementSpace");
504 
505  Array<int> count(flux.Size());
506  SumFluxAndCount(blfi, flux, count, wcoef, subdomain);
507 
508  if (ffes->Conforming()) // FIXME: nonconforming
509  {
510  // Accumulate flux and counts in parallel
511 
512  ffes->GroupComm().Reduce<double>(flux, GroupCommunicator::Sum);
513  ffes->GroupComm().Bcast<double>(flux);
514 
515  ffes->GroupComm().Reduce<int>(count, GroupCommunicator::Sum);
516  ffes->GroupComm().Bcast<int>(count);
517  }
518  else
519  {
520  MFEM_ABORT("Averaging on processor boundaries not implemented for "
521  "NC meshes yet.\n"
522  "Use L2ZZErrorEstimator() instead of ZZErrorEstimator().");
523  }
524 
525  // complete averaging
526  for (int i = 0; i < count.Size(); i++)
527  {
528  if (count[i] != 0) { flux(i) /= count[i]; }
529  }
530 
531  if (ffes->Nonconforming())
532  {
533  // On a partially conforming flux space, project on the conforming space.
534  // Using this code may lead to worse refinements in ex6, so we do not use
535  // it by default.
536 
537  // Vector conf_flux;
538  // flux.ConformingProject(conf_flux);
539  // flux.ConformingProlongate(conf_flux);
540  }
541 }
542 
543 
545  const ParGridFunction &x,
546  ParFiniteElementSpace &smooth_flux_fes,
547  ParFiniteElementSpace &flux_fes,
548  Vector &errors,
549  int norm_p, double solver_tol, int solver_max_it)
550 {
551  // Compute fluxes in discontinuous space
552  GridFunction flux(&flux_fes);
553  flux = 0.0;
554 
555  ParFiniteElementSpace *xfes = x.ParFESpace();
556  Array<int> xdofs, fdofs;
557  Vector el_x, el_f;
558 
559  for (int i = 0; i < xfes->GetNE(); i++)
560  {
561  xfes->GetElementVDofs(i, xdofs);
562  x.GetSubVector(xdofs, el_x);
563 
565  flux_integrator.ComputeElementFlux(*xfes->GetFE(i), *Transf, el_x,
566  *flux_fes.GetFE(i), el_f, false);
567 
568  flux_fes.GetElementVDofs(i, fdofs);
569  flux.AddElementVector(fdofs, el_f);
570  }
571 
572  // Assemble the linear system for L2 projection into the "smooth" space
573  ParBilinearForm *a = new ParBilinearForm(&smooth_flux_fes);
574  ParLinearForm *b = new ParLinearForm(&smooth_flux_fes);
576 
577  if (smooth_flux_fes.GetFE(0)->GetRangeType() == FiniteElement::SCALAR)
578  {
581  }
582  else
583  {
586  }
587 
588  b->Assemble();
589  a->Assemble();
590  a->Finalize();
591 
592  // The destination of the projected discontinuous flux
593  ParGridFunction smooth_flux(&smooth_flux_fes);
594  smooth_flux = 0.0;
595 
598  HypreParVector* X = smooth_flux.ParallelProject();
599 
600  delete a;
601  delete b;
602 
603  // Define and apply a parallel PCG solver for AX=B with the BoomerAMG
604  // preconditioner from hypre.
605  HypreBoomerAMG *amg = new HypreBoomerAMG(*A);
606  amg->SetPrintLevel(0);
607  HyprePCG *pcg = new HyprePCG(*A);
608  pcg->SetTol(solver_tol);
609  pcg->SetMaxIter(solver_max_it);
610  pcg->SetPrintLevel(0);
611  pcg->SetPreconditioner(*amg);
612  pcg->Mult(*B, *X);
613 
614  // Extract the parallel grid function corresponding to the finite element
615  // approximation X. This is the local solution on each processor.
616  smooth_flux = *X;
617 
618  delete A;
619  delete B;
620  delete X;
621  delete amg;
622  delete pcg;
623 
624  // Proceed through the elements one by one, and find the Lp norm differences
625  // between the flux as computed per element and the flux projected onto the
626  // smooth_flux_fes space.
627  double total_error = 0.0;
628  errors.SetSize(xfes->GetNE());
629  for (int i = 0; i < xfes->GetNE(); i++)
630  {
631  errors(i) = ComputeElementLpDistance(norm_p, i, smooth_flux, flux);
632  total_error += pow(errors(i), norm_p);
633  }
634 
635  double glob_error;
636  MPI_Allreduce(&total_error, &glob_error, 1, MPI_DOUBLE, MPI_SUM,
637  xfes->GetComm());
638 
639  return pow(glob_error, 1.0/norm_p);
640 }
641 
642 }
643 
644 #endif // MFEM_USE_MPI
int GetNFaceNeighbors() const
Definition: pmesh.hpp:154
void SetTol(double tol)
Definition: hypre.cpp:1993
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:175
int Size() const
Logical size of the array.
Definition: array.hpp:109
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:498
int GetVSize() const
Definition: fespace.hpp:163
int * GetJ()
Definition: table.hpp:108
void SetSpace(ParFiniteElementSpace *f)
Definition: pgridfunc.cpp:71
HypreParVector * NewTrueDofVector()
Definition: pfespace.hpp:203
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void SaveAsOne(std::ostream &out=std::cout)
Merge the local grid functions.
Definition: pgridfunc.cpp:329
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:155
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1169
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:138
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1013
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:310
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:133
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1361
ParFiniteElementSpace * pfes
Definition: pgridfunc.hpp:34
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:42
void Copy(Array &copy) const
Create a copy of the current array.
Definition: array.hpp:160
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Definition: vector.cpp:462
Delta function coefficient.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:106
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:584
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:213
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:314
void DivideByGroupSize(double *vec)
Scale a vector of true dofs.
Definition: pfespace.cpp:455
void ProjectDiscCoefficient(VectorCoefficient &coeff)
Definition: pgridfunc.cpp:270
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual const SparseMatrix * GetRestrictionMatrix()
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:244
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:249
double * GetData() const
Definition: vector.hpp:114
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:943
int Size_of_connections() const
Definition: table.hpp:92
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2008
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:842
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2185
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:184
void MakeRef(ParFiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new ParFiniteElementSpace.
Definition: pgridfunc.cpp:78
The BoomerAMG solver in hypre.
Definition: hypre.hpp:770
Class for parallel linear form.
Definition: plinearform.hpp:26
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:197
const FiniteElement * GetFaceNbrFE(int i) const
Definition: pfespace.cpp:876
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:136
void SetPrintLevel(int print_level)
Definition: hypre.hpp:811
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...
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:1035
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:68
void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:179
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:153
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:293
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add (element) subvector to the vector.
Definition: vector.cpp:550
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1998
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:72
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i&#39;th element.
Definition: fespace.hpp:205
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:461
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:630
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:475
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:210
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:544
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:92
void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:144
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:2515
bool Nonconforming() const
Definition: pfespace.hpp:261
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:117
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void Reduce(T *ldata, void(*Op)(OpData< T >))
Reduce within each group where the master is the root.
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:121
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2013
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void Destroy()
Destroy a vector.
Definition: vector.hpp:329
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1134
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:102
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:97
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1232
Class for parallel bilinear form.
void SetSpace(FiniteElementSpace *f)
Definition: gridfunc.cpp:171
int GetFaceNbrVSize() const
Definition: pfespace.hpp:249
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:176
Vector data type.
Definition: vector.hpp:36
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:177
Vector coefficient defined by a vector GridFunction.
int * GetI()
Definition: table.hpp:107
Class for parallel grid function.
Definition: pgridfunc.hpp:31
double * data
Definition: vector.hpp:41
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:174
HypreParMatrix * Dof_TrueDof_Matrix()
The true dof-to-dof interpolation matrix.
Definition: pfespace.hpp:193
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2034
int GetRangeType() const
Definition: fe.hpp:129
Class for parallel meshes.
Definition: pmesh.hpp:28
int GetFaceNbrRank(int fn) const
Definition: pmesh.cpp:1815
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.
void Bcast(T *data, int layout)
Broadcast within each group where the master is the root. The data layout can be: ...
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: pgridfunc.cpp:497
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:128
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:75
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const
Definition: fespace.cpp:68