MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pgridfunc.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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  /// Construct a ParGridFunction using a GridFunction as external data.
69  /** The parallel space @a *pf and the space used by @a *gf should match. The
70  data from @a *gf is used as the local data of the ParGridFunction on each
71  processor. The ParGridFunction does not assume ownership of the data. */
73 
74  /** @brief Creates grid function on (all) dofs from a given vector on the
75  true dofs, i.e. P tv. */
77 
78  /** @brief Construct a local ParGridFunction from the given *global*
79  GridFunction. If @a partitioning is NULL (default), the data from @a gf
80  is NOT copied. */
81  ParGridFunction(ParMesh *pmesh, const GridFunction *gf,
82  const int *partitioning = NULL);
83 
84  /** @brief Construct a ParGridFunction on a given ParMesh, @a pmesh, reading
85  from an std::istream.
86 
87  In the process, a ParFiniteElementSpace and a FiniteElementCollection are
88  constructed. The new ParGridFunction assumes ownership of both. */
89  ParGridFunction(ParMesh *pmesh, std::istream &input);
90 
91  /// Copy assignment. Only the data of the base class Vector is copied.
92  /** It is assumed that this object and @a rhs use ParFiniteElementSpace%s
93  that have the same size.
94 
95  @note Defining this method overwrites the implicitly defined copy
96  assignment operator. */
98  { return operator=((const Vector &)rhs); }
99 
100  /// Assign constant values to the ParGridFunction data.
101  ParGridFunction &operator=(double value)
102  { GridFunction::operator=(value); return *this; }
103 
104  /// Copy the data from a Vector to the ParGridFunction data.
106  { GridFunction::operator=(v); return *this; }
107 
108  ParFiniteElementSpace *ParFESpace() const { return pfes; }
109 
110  virtual void Update();
111 
112  /// Associate a new FiniteElementSpace with the ParGridFunction.
113  /** The ParGridFunction is resized using the SetSize() method. The new space
114  @a f is expected to be a ParFiniteElementSpace. */
115  virtual void SetSpace(FiniteElementSpace *f);
116 
117  /// Associate a new parallel space with the ParGridFunction.
119 
120  using GridFunction::MakeRef;
121 
122  /** @brief Make the ParGridFunction reference external data on a new
123  FiniteElementSpace. */
124  /** This method changes the FiniteElementSpace associated with the
125  ParGridFunction and sets the pointer @a v as external data in the
126  ParGridFunction. The new space @a f is expected to be a
127  ParFiniteElementSpace. */
128  virtual void MakeRef(FiniteElementSpace *f, double *v);
129 
130  /** @brief Make the ParGridFunction reference external data on a new
131  ParFiniteElementSpace. */
132  /** This method changes the ParFiniteElementSpace associated with the
133  ParGridFunction and sets the pointer @a v as external data in the
134  ParGridFunction. */
135  void MakeRef(ParFiniteElementSpace *f, double *v);
136 
137  /** @brief Make the ParGridFunction reference external data on a new
138  FiniteElementSpace. */
139  /** This method changes the FiniteElementSpace associated with the
140  ParGridFunction and sets the data of the Vector @a v (plus the @a
141  v_offset) as external data in the ParGridFunction. The new space @a f is
142  expected to be a ParFiniteElementSpace.
143  @note This version of the method will also perform bounds checks when
144  the build option MFEM_DEBUG is enabled. */
145  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
146 
147  /** @brief Make the ParGridFunction reference external data on a new
148  ParFiniteElementSpace. */
149  /** This method changes the ParFiniteElementSpace associated with the
150  ParGridFunction and sets the data of the Vector @a v (plus the
151  @a v_offset) as external data in the ParGridFunction.
152  @note This version of the method will also perform bounds checks when
153  the build option MFEM_DEBUG is enabled. */
154  void MakeRef(ParFiniteElementSpace *f, Vector &v, int v_offset);
155 
156  /** Set the grid function on (all) dofs from a given vector on the
157  true dofs, i.e. P tv. */
158  void Distribute(const Vector *tv);
159  void Distribute(const Vector &tv) { Distribute(&tv); }
160  void AddDistribute(double a, const Vector *tv);
161  void AddDistribute(double a, const Vector &tv) { AddDistribute(a, &tv); }
162 
163  /// Set the GridFunction from the given true-dof vector.
164  virtual void SetFromTrueDofs(const Vector &tv) { Distribute(tv); }
165 
166  /// Short semantic for Distribute()
168  { Distribute(&tv); return (*this); }
169 
171 
172  /// Returns the true dofs in a new HypreParVector
173  HypreParVector *GetTrueDofs() const;
174 
175  /// Returns the vector averaged on the true dofs.
176  void ParallelAverage(Vector &tv) const;
177 
178  /// Returns the vector averaged on the true dofs.
179  void ParallelAverage(HypreParVector &tv) const;
180 
181  /// Returns a new vector averaged on the true dofs.
183 
184  /// Returns the vector restricted to the true dofs.
185  void ParallelProject(Vector &tv) const;
186 
187  /// Returns the vector restricted to the true dofs.
188  void ParallelProject(HypreParVector &tv) const;
189 
190  /// Returns a new vector restricted to the true dofs.
192 
193  /// Returns the vector assembled on the true dofs.
194  void ParallelAssemble(Vector &tv) const;
195 
196  /// Returns the vector assembled on the true dofs.
197  void ParallelAssemble(HypreParVector &tv) const;
198 
199  /// Returns a new vector assembled on the true dofs.
201 
202  void ExchangeFaceNbrData();
204  const Vector &FaceNbrData() const { return face_nbr_data; }
205 
206  // Redefine to handle the case when i is a face-neighbor element
207  virtual double GetValue(int i, const IntegrationPoint &ip,
208  int vdim = 1) const;
210  { return GetValue(T.ElementNo, T.GetIntPoint()); }
211 
212  // Redefine to handle the case when T describes a face-neighbor element
213  virtual double GetValue(ElementTransformation &T, const IntegrationPoint &ip,
214  int comp = 0, Vector *tr = NULL) const;
215 
216  virtual void GetVectorValue(int i, const IntegrationPoint &ip,
217  Vector &val) const;
218 
219  // Redefine to handle the case when T describes a face-neighbor element
220  virtual void GetVectorValue(ElementTransformation &T,
221  const IntegrationPoint &ip,
222  Vector &val, Vector *tr = NULL) const;
223 
224  /** Sets the output vector @a dof_vals to the values of the degrees of
225  freedom of element @a el. If @a el is greater than or equal to the number
226  of local elements, it will be interpreted as a shifted index of a face
227  neighbor element. */
228  virtual void GetElementDofValues(int el, Vector &dof_vals) const;
229 
231  virtual void ProjectCoefficient(Coefficient &coeff);
232 
234  /** @brief Project a discontinuous vector coefficient as a grid function on
235  a continuous finite element space. The values in shared dofs are
236  determined from the element with maximal attribute. */
237  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
238 
239  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
240 
241  virtual void ProjectDiscCoefficient(VectorCoefficient &vcoeff, AvgType type);
242 
244 
245  // Only the values in the master are guaranteed to be correct!
247  Array<int> &attr)
248  { ProjectBdrCoefficient(NULL, &vcoeff, attr); }
249 
250  // Only the values in the master are guaranteed to be correct!
251  virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr)
252  { ProjectBdrCoefficient(coeff, NULL, attr); }
253 
254  // Only the values in the master are guaranteed to be correct!
255  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
256  Array<int> &bdr_attr);
257 
258  virtual double ComputeL1Error(Coefficient *exsol[],
259  const IntegrationRule *irs[] = NULL) const
260  {
262  *exsol, NULL, 1, NULL, irs), pfes->GetComm());
263  }
264 
265  virtual double ComputeL1Error(Coefficient &exsol,
266  const IntegrationRule *irs[] = NULL) const
267  { return ComputeLpError(1.0, exsol, NULL, irs); }
268 
269  virtual double ComputeL1Error(VectorCoefficient &exsol,
270  const IntegrationRule *irs[] = NULL) const
271  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
272 
273  virtual double ComputeL2Error(Coefficient *exsol[],
274  const IntegrationRule *irs[] = NULL) const
275  {
276  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs),
277  pfes->GetComm());
278  }
279 
280  virtual double ComputeL2Error(Coefficient &exsol,
281  const IntegrationRule *irs[] = NULL) const
282  { return ComputeLpError(2.0, exsol, NULL, irs); }
283 
284  virtual double ComputeL2Error(VectorCoefficient &exsol,
285  const IntegrationRule *irs[] = NULL,
286  Array<int> *elems = NULL) const
287  {
288  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
289  pfes->GetComm());
290  }
291 
292  /// Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements
293  virtual double ComputeGradError(VectorCoefficient *exgrad,
294  const IntegrationRule *irs[] = NULL) const
295  {
296  return GlobalLpNorm(2.0, GridFunction::ComputeGradError(exgrad,irs),
297  pfes->GetComm());
298  }
299 
300  /// Returns ||curl u_ex - curl u_h||_L2 for ND elements
301  virtual double ComputeCurlError(VectorCoefficient *excurl,
302  const IntegrationRule *irs[] = NULL) const
303  {
304  return GlobalLpNorm(2.0, GridFunction::ComputeCurlError(excurl,irs),
305  pfes->GetComm());
306  }
307 
308  /// Returns ||div u_ex - div u_h||_L2 for RT elements
309  virtual double ComputeDivError(Coefficient *exdiv,
310  const IntegrationRule *irs[] = NULL) const
311  {
312  return GlobalLpNorm(2.0, GridFunction::ComputeDivError(exdiv,irs),
313  pfes->GetComm());
314  }
315 
316  /// Returns the Face Jumps error for L2 elements
317  virtual double ComputeDGFaceJumpError(Coefficient *exsol,
318  Coefficient *ell_coeff,
319  JumpScaling jump_scaling,
320  const IntegrationRule *irs[]=NULL)
321  const;
322 
323  /// Returns either the H1-seminorm or the DG Face Jumps error or both
324  /// depending on norm_type = 1, 2, 3
325  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
326  Coefficient *ell_coef, double Nu,
327  int norm_type) const
328  {
329  return GlobalLpNorm(2.0,
330  GridFunction::ComputeH1Error(exsol,exgrad,ell_coef,
331  Nu, norm_type),
332  pfes->GetComm());
333  }
334 
335  /// Returns the error measured in H1-norm for H1 elements or in "broken"
336  /// H1-norm for L2 elements
337  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
338  const IntegrationRule *irs[] = NULL) const
339  {
340  return GlobalLpNorm(2.0, GridFunction::ComputeH1Error(exsol,exgrad,irs),
341  pfes->GetComm());
342  }
343 
344  /// Returns the error measured H(div)-norm for RT elements
345  virtual double ComputeHDivError(VectorCoefficient *exsol,
346  Coefficient *exdiv,
347  const IntegrationRule *irs[] = NULL) const
348  {
349  return GlobalLpNorm(2.0, GridFunction::ComputeHDivError(exsol,exdiv,irs),
350  pfes->GetComm());
351  }
352 
353  /// Returns the error measured H(curl)-norm for ND elements
354  virtual double ComputeHCurlError(VectorCoefficient *exsol,
355  VectorCoefficient *excurl,
356  const IntegrationRule *irs[] = NULL) const
357  {
358  return GlobalLpNorm(2.0,
359  GridFunction::ComputeHCurlError(exsol,excurl,irs),
360  pfes->GetComm());
361  }
362 
363  virtual double ComputeMaxError(Coefficient *exsol[],
364  const IntegrationRule *irs[] = NULL) const
365  {
366  return GlobalLpNorm(infinity(),
367  GridFunction::ComputeMaxError(exsol, irs),
368  pfes->GetComm());
369  }
370 
371  virtual double ComputeMaxError(Coefficient &exsol,
372  const IntegrationRule *irs[] = NULL) const
373  {
374  return ComputeLpError(infinity(), exsol, NULL, irs);
375  }
376 
377  virtual double ComputeMaxError(VectorCoefficient &exsol,
378  const IntegrationRule *irs[] = NULL) const
379  {
380  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
381  }
382 
383  virtual double ComputeLpError(const double p, Coefficient &exsol,
384  Coefficient *weight = NULL,
385  const IntegrationRule *irs[] = NULL) const
386  {
388  p, exsol, weight, irs), pfes->GetComm());
389  }
390 
391  /** When given a vector weight, compute the pointwise (scalar) error as the
392  dot product of the vector error with the vector weight. Otherwise, the
393  scalar error is the l_2 norm of the vector error. */
394  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
395  Coefficient *weight = NULL,
396  VectorCoefficient *v_weight = NULL,
397  const IntegrationRule *irs[] = NULL) const
398  {
400  p, exsol, weight, v_weight, irs), pfes->GetComm());
401  }
402 
403  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
404  GridFunction &flux,
405  bool wcoef = true, int subdomain = -1);
406 
407  /** Save the local portion of the ParGridFunction. This differs from the
408  serial GridFunction::Save in that it takes into account the signs of
409  the local dofs. */
410  virtual void Save(std::ostream &out) const;
411 
412  /// Save the ParGridFunction to a single file (written using MPI rank 0). The
413  /// given @a precision will be used for ASCII output.
414  void SaveAsOne(const char *fname, int precision=16) const;
415 
416  /// Save the ParGridFunction to files (one for each MPI rank). The files will
417  /// be given suffixes according to the MPI rank. The given @a precision will
418  /// be used for ASCII output.
419  virtual void Save(const char *fname, int precision=16) const;
420 
421 #ifdef MFEM_USE_ADIOS2
422  /** Save the local portion of the ParGridFunction. This differs from the
423  serial GridFunction::Save in that it takes into account the signs of
424  the local dofs. */
425  virtual void Save(
426  adios2stream &out, const std::string &variable_name,
428 #endif
429 
430  /// Merge the local grid functions
431  void SaveAsOne(std::ostream &out = mfem::out) const;
432 
433  virtual ~ParGridFunction() { }
434 };
435 
436 
437 /** Performs a global L2 projection (through a HypreBoomerAMG solve) of flux
438  from supplied discontinuous space into supplied smooth (continuous, or at
439  least conforming) space, and computes the Lp norms of the differences
440  between them on each element. This is one approach to handling conforming
441  and non-conforming elements in parallel. Returns the total error estimate. */
442 double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator,
443  const ParGridFunction &x,
444  ParFiniteElementSpace &smooth_flux_fes,
445  ParFiniteElementSpace &flux_fes,
446  Vector &errors, int norm_p = 2, double solver_tol = 1e-12,
447  int solver_max_it = 200);
448 
449 }
450 
451 #endif // MFEM_USE_MPI
452 
453 #endif
virtual double ComputeLpError(const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:394
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Memory< double > data
Definition: vector.hpp:64
void AddDistribute(double a, const Vector &tv)
Definition: pgridfunc.hpp:161
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
Definition: pgridfunc.hpp:53
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:198
Base class for vector Coefficients that optionally depend on time and space.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:456
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:181
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: gridfunc.cpp:2749
virtual void GetElementDofValues(int el, Vector &dof_vals) const
Definition: pgridfunc.cpp:474
virtual ~ParGridFunction()
Definition: pgridfunc.hpp:433
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2433
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:167
ParGridFunction(ParFiniteElementSpace *pf)
Definition: pgridfunc.hpp:56
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:264
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:841
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Definition: pgridfunc.cpp:514
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:87
Abstract parallel finite element space.
Definition: pfespace.hpp:28
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:2976
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:493
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
Definition: pgridfunc.hpp:44
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:90
void Distribute(const Vector &tv)
Definition: pgridfunc.hpp:159
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:246
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: pgridfunc.hpp:301
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:345
double GetValue(ElementTransformation &T)
Definition: pgridfunc.hpp:209
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
MPI_Comm GetComm() const
Definition: pfespace.hpp:263
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
ParGridFunction & operator=(double value)
Assign constant values to the ParGridFunction data.
Definition: pgridfunc.hpp:101
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:164
const Vector & FaceNbrData() const
Definition: pgridfunc.hpp:204
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:2950
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:269
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:677
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:201
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition: gridfunc.cpp:332
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:383
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: gridfunc.hpp:112
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3156
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
Definition: pgridfunc.hpp:105
Vector & FaceNbrData()
Definition: pgridfunc.hpp:203
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: pgridfunc.cpp:321
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:99
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:258
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:2709
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:1033
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: pgridfunc.hpp:97
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:280
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:601
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: pgridfunc.cpp:1068
virtual double ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:363
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:520
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:1105
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:132
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:293
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:867
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: gridfunc.cpp:2791
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:371
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:377
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:2985
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:251
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:25
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:273
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:337
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:639
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:164
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:138
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2278
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:265
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3050
virtual double ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:284
Vector data type.
Definition: vector.hpp:60
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: pgridfunc.hpp:325
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: pgridfunc.hpp:309
Class for parallel grid function.
Definition: pgridfunc.hpp:32
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Class for parallel meshes.
Definition: pmesh.hpp:32
ParGridFunction(ParFiniteElementSpace *pf, double *data)
Construct a ParGridFunction using previously allocated array data.
Definition: pgridfunc.hpp:65
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:354
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:424
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108