MFEM  v4.6.0
Finite element discretization library
pgridfunc.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 #ifndef MFEM_PGRIDFUNC
13 #define MFEM_PGRIDFUNC
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include "../general/globals.hpp"
20 #include "pfespace.hpp"
21 #include "gridfunc.hpp"
22 #include <iostream>
23 #include <limits>
24 
25 namespace mfem
26 {
27 
28 /// Compute a global Lp norm from the local Lp norms computed by each processor
29 double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm);
30 
31 /// Class for parallel grid function
33 {
34 protected:
35  ParFiniteElementSpace *pfes; ///< Points to the same object as #fes
36 
37  /** @brief Vector used to store data from face-neighbor processors,
38  initialized by ExchangeFaceNbrData(). */
40 
41  /** @brief Vector used as an MPI buffer to send face-neighbor data
42  in ExchangeFaceNbrData() to neighboring processors. */
43  //TODO: Use temporary memory to avoid CUDA malloc allocation cost.
45 
47  Array<int> &attr);
48 
49 public:
50  ParGridFunction() { pfes = NULL; }
51 
52  /// Copy constructor. The internal vector #face_nbr_data is not copied.
54  : GridFunction(orig), pfes(orig.pfes) { }
55 
57 
58  /// Construct a ParGridFunction using previously allocated array @a data.
59  /** The ParGridFunction does not assume ownership of @a data which is assumed
60  to be of size at least `pf->GetVSize()`. Similar to the GridFunction and
61  Vector constructors for externally allocated array, the pointer @a data
62  can be NULL. The data array can be replaced later using the method
63  SetData().
64  */
66  GridFunction(pf, data), pfes(pf) { }
67 
68  /** @brief Construct a ParGridFunction using previously allocated Vector
69  @a base starting at the given offset, @a base_offset. */
70  ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset = 0)
71  : GridFunction(pf, base, base_offset), pfes(pf) { }
72 
73  /// Construct a ParGridFunction using a GridFunction as external data.
74  /** The parallel space @a *pf and the space used by @a *gf should match. The
75  data from @a *gf is used as the local data of the ParGridFunction on each
76  processor. The ParGridFunction does not assume ownership of the data. */
78 
79  /** @brief Creates grid function on (all) dofs from a given vector on the
80  true dofs, i.e. P tv. */
82 
83  /** @brief Construct a local ParGridFunction from the given *global*
84  GridFunction. If @a partitioning is NULL (default), the data from @a gf
85  is NOT copied. */
86  ParGridFunction(ParMesh *pmesh, const GridFunction *gf,
87  const int *partitioning = NULL);
88 
89  /** @brief Construct a ParGridFunction on a given ParMesh, @a pmesh, reading
90  from an std::istream.
91 
92  In the process, a ParFiniteElementSpace and a FiniteElementCollection are
93  constructed. The new ParGridFunction assumes ownership of both. */
94  ParGridFunction(ParMesh *pmesh, std::istream &input);
95 
96  /// Copy assignment. Only the data of the base class Vector is copied.
97  /** It is assumed that this object and @a rhs use ParFiniteElementSpace%s
98  that have the same size.
99 
100  @note Defining this method overwrites the implicitly defined copy
101  assignment operator. */
103  { return operator=((const Vector &)rhs); }
104 
105  /// Assign constant values to the ParGridFunction data.
106  ParGridFunction &operator=(double value)
107  { GridFunction::operator=(value); return *this; }
108 
109  /// Copy the data from a Vector to the ParGridFunction data.
111  { GridFunction::operator=(v); return *this; }
112 
113  ParFiniteElementSpace *ParFESpace() const { return pfes; }
114 
115  virtual void Update();
116 
117  /// Associate a new FiniteElementSpace with the ParGridFunction.
118  /** The ParGridFunction is resized using the SetSize() method. The new space
119  @a f is expected to be a ParFiniteElementSpace. */
120  virtual void SetSpace(FiniteElementSpace *f);
121 
122  /// Associate a new parallel space with the ParGridFunction.
124 
125  using GridFunction::MakeRef;
126 
127  /** @brief Make the ParGridFunction reference external data on a new
128  FiniteElementSpace. */
129  /** This method changes the FiniteElementSpace associated with the
130  ParGridFunction and sets the pointer @a v as external data in the
131  ParGridFunction. The new space @a f is expected to be a
132  ParFiniteElementSpace. */
133  virtual void MakeRef(FiniteElementSpace *f, double *v);
134 
135  /** @brief Make the ParGridFunction reference external data on a new
136  ParFiniteElementSpace. */
137  /** This method changes the ParFiniteElementSpace associated with the
138  ParGridFunction and sets the pointer @a v as external data in the
139  ParGridFunction. */
140  void MakeRef(ParFiniteElementSpace *f, double *v);
141 
142  /** @brief Make the ParGridFunction reference external data on a new
143  FiniteElementSpace. */
144  /** This method changes the FiniteElementSpace associated with the
145  ParGridFunction and sets the data of the Vector @a v (plus the @a
146  v_offset) as external data in the ParGridFunction. The new space @a f is
147  expected to be a ParFiniteElementSpace.
148  @note This version of the method will also perform bounds checks when
149  the build option MFEM_DEBUG is enabled. */
150  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
151 
152  /** @brief Make the ParGridFunction reference external data on a new
153  ParFiniteElementSpace. */
154  /** This method changes the ParFiniteElementSpace associated with the
155  ParGridFunction and sets the data of the Vector @a v (plus the
156  @a v_offset) as external data in the ParGridFunction.
157  @note This version of the method will also perform bounds checks when
158  the build option MFEM_DEBUG is enabled. */
159  void MakeRef(ParFiniteElementSpace *f, Vector &v, int v_offset);
160 
161  /** Set the grid function on (all) dofs from a given vector on the
162  true dofs, i.e. P tv. */
163  void Distribute(const Vector *tv);
164  void Distribute(const Vector &tv) { Distribute(&tv); }
165  void AddDistribute(double a, const Vector *tv);
166  void AddDistribute(double a, const Vector &tv) { AddDistribute(a, &tv); }
167 
168  /// Set the GridFunction from the given true-dof vector.
169  virtual void SetFromTrueDofs(const Vector &tv) { Distribute(tv); }
170 
171  /// Short semantic for Distribute()
173  { Distribute(&tv); return (*this); }
174 
176 
177  /// Returns the true dofs in a new HypreParVector
178  HypreParVector *GetTrueDofs() const;
179 
180  /// Returns the vector averaged on the true dofs.
181  void ParallelAverage(Vector &tv) const;
182 
183  /// Returns the vector averaged on the true dofs.
184  void ParallelAverage(HypreParVector &tv) const;
185 
186  /// Returns a new vector averaged on the true dofs.
188 
189  /// Returns the vector restricted to the true dofs.
190  void ParallelProject(Vector &tv) const;
191 
192  /// Returns the vector restricted to the true dofs.
193  void ParallelProject(HypreParVector &tv) const;
194 
195  /// Returns a new vector restricted to the true dofs.
197 
198  /// Returns the vector assembled on the true dofs.
199  void ParallelAssemble(Vector &tv) const;
200 
201  /// Returns the vector assembled on the true dofs.
202  void ParallelAssemble(HypreParVector &tv) const;
203 
204  /// Returns a new vector assembled on the true dofs.
206 
207  void ExchangeFaceNbrData();
209  const Vector &FaceNbrData() const { return face_nbr_data; }
210 
211  // Redefine to handle the case when i is a face-neighbor element
212  virtual double GetValue(int i, const IntegrationPoint &ip,
213  int vdim = 1) const;
215  { return GetValue(T.ElementNo, T.GetIntPoint()); }
216 
217  // Redefine to handle the case when T describes a face-neighbor element
218  virtual double GetValue(ElementTransformation &T, const IntegrationPoint &ip,
219  int comp = 0, Vector *tr = NULL) const;
220 
221  virtual void GetVectorValue(int i, const IntegrationPoint &ip,
222  Vector &val) const;
223 
224  // Redefine to handle the case when T describes a face-neighbor element
225  virtual void GetVectorValue(ElementTransformation &T,
226  const IntegrationPoint &ip,
227  Vector &val, Vector *tr = NULL) const;
228 
229  /** @brief For each vdof, counts how many elements contain the vdof,
230  as containment is determined by FiniteElementSpace::GetElementVDofs(). */
231  virtual void CountElementsPerVDof(Array<int> &elem_per_vdof) const;
232 
233  /// Parallel version of GridFunction::GetDerivative(); see its documentation.
234  void GetDerivative(int comp, int der_comp, ParGridFunction &der);
235 
236  /** Sets the output vector @a dof_vals to the values of the degrees of
237  freedom of element @a el. If @a el is greater than or equal to the number
238  of local elements, it will be interpreted as a shifted index of a face
239  neighbor element. */
240  virtual void GetElementDofValues(int el, Vector &dof_vals) const;
241 
243  virtual void ProjectCoefficient(Coefficient &coeff);
244 
246  /** @brief Project a discontinuous vector coefficient as a grid function on
247  a continuous finite element space. The values in shared dofs are
248  determined from the element with maximal attribute. */
249  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
250 
251  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
252 
253  virtual void ProjectDiscCoefficient(VectorCoefficient &vcoeff, AvgType type);
254 
256 
257  // Only the values in the master are guaranteed to be correct!
259  Array<int> &attr)
260  { ProjectBdrCoefficient(NULL, &vcoeff, attr); }
261 
262  // Only the values in the master are guaranteed to be correct!
263  virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr)
264  { ProjectBdrCoefficient(coeff, NULL, attr); }
265 
266  // Only the values in the master are guaranteed to be correct!
267  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
268  Array<int> &bdr_attr);
269 
270  virtual double ComputeL1Error(Coefficient *exsol[],
271  const IntegrationRule *irs[] = NULL) const
272  {
274  *exsol, NULL, 1, NULL, irs), pfes->GetComm());
275  }
276 
277  virtual double ComputeL1Error(Coefficient &exsol,
278  const IntegrationRule *irs[] = NULL) const
279  { return ComputeLpError(1.0, exsol, NULL, irs); }
280 
281  virtual double ComputeL1Error(VectorCoefficient &exsol,
282  const IntegrationRule *irs[] = NULL) const
283  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
284 
285  virtual double ComputeL2Error(Coefficient *exsol[],
286  const IntegrationRule *irs[] = NULL,
287  const Array<int> *elems = NULL) const
288  {
289  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
290  pfes->GetComm());
291  }
292 
293  virtual double ComputeL2Error(Coefficient &exsol,
294  const IntegrationRule *irs[] = NULL,
295  const Array<int> *elems = NULL) const
296  {
297  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
298  pfes->GetComm());
299  }
300 
301 
302  virtual double ComputeL2Error(VectorCoefficient &exsol,
303  const IntegrationRule *irs[] = NULL,
304  const Array<int> *elems = NULL) const
305  {
306  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
307  pfes->GetComm());
308  }
309 
310  /// Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements
311  virtual double ComputeGradError(VectorCoefficient *exgrad,
312  const IntegrationRule *irs[] = NULL) const
313  {
314  return GlobalLpNorm(2.0, GridFunction::ComputeGradError(exgrad,irs),
315  pfes->GetComm());
316  }
317 
318  /// Returns ||curl u_ex - curl u_h||_L2 for ND elements
319  virtual double ComputeCurlError(VectorCoefficient *excurl,
320  const IntegrationRule *irs[] = NULL) const
321  {
322  return GlobalLpNorm(2.0, GridFunction::ComputeCurlError(excurl,irs),
323  pfes->GetComm());
324  }
325 
326  /// Returns ||div u_ex - div u_h||_L2 for RT elements
327  virtual double ComputeDivError(Coefficient *exdiv,
328  const IntegrationRule *irs[] = NULL) const
329  {
330  return GlobalLpNorm(2.0, GridFunction::ComputeDivError(exdiv,irs),
331  pfes->GetComm());
332  }
333 
334  /// Returns the Face Jumps error for L2 elements
335  virtual double ComputeDGFaceJumpError(Coefficient *exsol,
336  Coefficient *ell_coeff,
337  JumpScaling jump_scaling,
338  const IntegrationRule *irs[]=NULL)
339  const;
340 
341  /// Returns either the H1-seminorm or the DG Face Jumps error or both
342  /// depending on norm_type = 1, 2, 3
343  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
344  Coefficient *ell_coef, double Nu,
345  int norm_type) const
346  {
347  return GlobalLpNorm(2.0,
348  GridFunction::ComputeH1Error(exsol,exgrad,ell_coef,
349  Nu, norm_type),
350  pfes->GetComm());
351  }
352 
353  /// Returns the error measured in H1-norm for H1 elements or in "broken"
354  /// H1-norm for L2 elements
355  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
356  const IntegrationRule *irs[] = NULL) const
357  {
358  return GlobalLpNorm(2.0, GridFunction::ComputeH1Error(exsol,exgrad,irs),
359  pfes->GetComm());
360  }
361 
362  /// Returns the error measured H(div)-norm for RT elements
363  virtual double ComputeHDivError(VectorCoefficient *exsol,
364  Coefficient *exdiv,
365  const IntegrationRule *irs[] = NULL) const
366  {
367  return GlobalLpNorm(2.0, GridFunction::ComputeHDivError(exsol,exdiv,irs),
368  pfes->GetComm());
369  }
370 
371  /// Returns the error measured H(curl)-norm for ND elements
372  virtual double ComputeHCurlError(VectorCoefficient *exsol,
373  VectorCoefficient *excurl,
374  const IntegrationRule *irs[] = NULL) const
375  {
376  return GlobalLpNorm(2.0,
377  GridFunction::ComputeHCurlError(exsol,excurl,irs),
378  pfes->GetComm());
379  }
380 
381  virtual double ComputeMaxError(Coefficient *exsol[],
382  const IntegrationRule *irs[] = NULL) const
383  {
384  return GlobalLpNorm(infinity(),
385  GridFunction::ComputeMaxError(exsol, irs),
386  pfes->GetComm());
387  }
388 
389  virtual double ComputeMaxError(Coefficient &exsol,
390  const IntegrationRule *irs[] = NULL) const
391  {
392  return ComputeLpError(infinity(), exsol, NULL, irs);
393  }
394 
395  virtual double ComputeMaxError(VectorCoefficient &exsol,
396  const IntegrationRule *irs[] = NULL) const
397  {
398  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
399  }
400 
401  virtual double ComputeLpError(const double p, Coefficient &exsol,
402  Coefficient *weight = NULL,
403  const IntegrationRule *irs[] = NULL,
404  const Array<int> *elems = NULL) const
405  {
406  return GlobalLpNorm(p, GridFunction::ComputeLpError(p, exsol, weight, irs,
407  elems), pfes->GetComm());
408  }
409 
410  /** When given a vector weight, compute the pointwise (scalar) error as the
411  dot product of the vector error with the vector weight. Otherwise, the
412  scalar error is the l_2 norm of the vector error. */
413  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
414  Coefficient *weight = NULL,
415  VectorCoefficient *v_weight = NULL,
416  const IntegrationRule *irs[] = NULL) const
417  {
419  p, exsol, weight, v_weight, irs), pfes->GetComm());
420  }
421 
422  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
423  GridFunction &flux,
424  bool wcoef = true, int subdomain = -1);
425 
426  /** Save the local portion of the ParGridFunction. This differs from the
427  serial GridFunction::Save in that it takes into account the signs of
428  the local dofs. */
429  virtual void Save(std::ostream &out) const;
430 
431  /// Save the ParGridFunction to a single file (written using MPI rank 0). The
432  /// given @a precision will be used for ASCII output.
433  void SaveAsOne(const char *fname, int precision=16) const;
434 
435  /// Save the ParGridFunction to files (one for each MPI rank). The files will
436  /// be given suffixes according to the MPI rank. The given @a precision will
437  /// be used for ASCII output.
438  virtual void Save(const char *fname, int precision=16) const;
439 
440  /// Returns a GridFunction on MPI rank @a save_rank that does not have any
441  /// duplication of vertices/nodes at processor boundaries.
442  /// serial_mesh is obtained using ParMesh::GetSerialMesh(save_rank).
443  /// Note that the @ save_rank argument must match for the
444  /// ParMesh::GetSerialMesh and GetSerialGridFunction method.
445  GridFunction GetSerialGridFunction(int save_rank, Mesh &serial_mesh) const;
446 
447  /// Write the serial GridFunction a single file (written using MPI rank 0).
448  /// The given @a precision will be used for ASCII output.
449  void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const;
450 
451 #ifdef MFEM_USE_ADIOS2
452  /** Save the local portion of the ParGridFunction. This differs from the
453  serial GridFunction::Save in that it takes into account the signs of
454  the local dofs. */
455  virtual void Save(
456  adios2stream &out, const std::string &variable_name,
458 #endif
459 
460  /// Merge the local grid functions
461  void SaveAsOne(std::ostream &out = mfem::out) const;
462 
463  virtual ~ParGridFunction() { }
464 };
465 
466 
467 /** Performs a global L2 projection (through a HypreBoomerAMG solve) of flux
468  from supplied discontinuous space into supplied smooth (continuous, or at
469  least conforming) space, and computes the Lp norms of the differences
470  between them on each element. This is one approach to handling conforming
471  and non-conforming elements in parallel. Returns the total error estimate. */
472 double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator,
473  const ParGridFunction &x,
474  ParFiniteElementSpace &smooth_flux_fes,
475  ParFiniteElementSpace &flux_fes,
476  Vector &errors, int norm_p = 2, double solver_tol = 1e-12,
477  int solver_max_it = 200);
478 
479 }
480 
481 #endif // MFEM_USE_MPI
482 
483 #endif
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:2786
ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset=0)
Construct a ParGridFunction using previously allocated Vector base starting at the given offset...
Definition: pgridfunc.hpp:70
virtual void GetElementDofValues(int el, Vector &dof_vals) const
Definition: pgridfunc.cpp:542
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
Definition: gridfunc.cpp:2920
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:355
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void GetDerivative(int comp, int der_comp, ParGridFunction &der)
Parallel version of GridFunction::GetDerivative(); see its documentation.
Definition: pgridfunc.cpp:521
Memory< double > data
Definition: vector.hpp:62
void AddDistribute(double a, const Vector &tv)
Definition: pgridfunc.hpp:166
virtual double ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:302
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
Definition: pgridfunc.hpp:53
Base class for vector Coefficients that optionally depend on time and space.
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:935
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:207
virtual ~ParGridFunction()
Definition: pgridfunc.hpp:463
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2589
ParFiniteElementSpace * pfes
Points to the same object as fes.
Definition: pgridfunc.hpp:35
ParGridFunction & operator=(const HypreParVector &tv)
Short semantic for Distribute()
Definition: pgridfunc.hpp:172
ParGridFunction(ParFiniteElementSpace *pf)
Definition: pgridfunc.hpp:56
virtual void CountElementsPerVDof(Array< int > &elem_per_vdof) const
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
Definition: pgridfunc.cpp:512
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Definition: pgridfunc.cpp:582
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:96
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
Definition: pgridfunc.hpp:44
void Distribute(const Vector &tv)
Definition: pgridfunc.hpp:164
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:270
virtual void ProjectBdrCoefficient(VectorCoefficient &vcoeff, Array< int > &attr)
Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with ...
Definition: pgridfunc.hpp:258
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:401
virtual double ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:381
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
double GetValue(ElementTransformation &T)
Definition: pgridfunc.hpp:214
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:111
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:574
const Vector & FaceNbrData() const
Definition: pgridfunc.hpp:209
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition: gridfunc.cpp:366
ParGridFunction & operator=(double value)
Assign constant values to the ParGridFunction data.
Definition: pgridfunc.hpp:106
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:169
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:90
virtual double ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
Definition: gridfunc.cpp:3186
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: gridfunc.cpp:2960
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: pgridfunc.hpp:327
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:190
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:203
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: pgridfunc.hpp:343
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured H(curl)-norm for ND elements.
Definition: pgridfunc.hpp:372
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: gridfunc.hpp:116
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3260
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:285
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
Definition: pgridfunc.hpp:110
Vector & FaceNbrData()
Definition: pgridfunc.hpp:208
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:161
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:1198
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:277
virtual double ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured H(div)-norm for RT elements.
Definition: pgridfunc.hpp:363
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: pgridfunc.hpp:102
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:669
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:173
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: pgridfunc.cpp:1233
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:3160
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
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: gridfunc.cpp:3001
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
Returns the Face Jumps error for L2 elements.
Definition: pgridfunc.cpp:745
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:389
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:1270
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:141
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:273
virtual double ComputeLpError(const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:413
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: pgridfunc.cpp:341
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:909
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
Definition: gridfunc.cpp:3195
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:3366
double a
Definition: lissajous.cpp:41
virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array< int > &attr)
Project a set of Coefficients on the components of the GridFunction, modifying only DOFs on the bound...
Definition: pgridfunc.hpp:263
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
Definition: pgridfunc.hpp:39
Class for integration point with weight.
Definition: intrules.hpp:31
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:281
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:707
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:147
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2415
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
Definition: pgridfunc.hpp:311
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||u_ex - u_h||_L2 for H1 or L2 elements.
Definition: pgridfunc.hpp:293
MPI_Comm GetComm() const
Definition: pfespace.hpp:269
void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const
Definition: pgridfunc.cpp:947
Vector data type.
Definition: vector.hpp:58
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:152
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: pgridfunc.hpp:319
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Class for parallel meshes.
Definition: pmesh.hpp:32
GridFunction GetSerialGridFunction(int save_rank, Mesh &serial_mesh) const
Definition: pgridfunc.cpp:961
ParGridFunction(ParFiniteElementSpace *pf, double *data)
Construct a ParGridFunction using previously allocated array data.
Definition: pgridfunc.hpp:65
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:468
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:395