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