MFEM  v4.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, 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 #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 
42  Array<int> &attr);
43 
44 public:
45  ParGridFunction() { pfes = NULL; }
46 
47  /// Copy constructor. The internal vector #face_nbr_data is not copied.
49  : GridFunction(orig), pfes(orig.pfes) { }
50 
52 
53  /// Construct a ParGridFunction using previously allocated array @a data.
54  /** The ParGridFunction does not assume ownership of @a data which is assumed
55  to be of size at least `pf->GetVSize()`. Similar to the GridFunction and
56  Vector constructors for externally allocated array, the pointer @a data
57  can be NULL. The data array can be replaced later using the method
58  SetData().
59  */
61  GridFunction(pf, data), pfes(pf) { }
62 
63  /// Construct a ParGridFunction using a GridFunction as external data.
64  /** The parallel space @a *pf and the space used by @a *gf should match. The
65  data from @a *gf is used as the local data of the ParGridFunction on each
66  processor. The ParGridFunction does not assume ownership of the data. */
68 
69  /** @brief Creates grid function on (all) dofs from a given vector on the
70  true dofs, i.e. P tv. */
72 
73  /** @brief Construct a local ParGridFunction from the given *global*
74  GridFunction. If @a partitioning is NULL (default), the data from @a gf
75  is NOT copied. */
76  ParGridFunction(ParMesh *pmesh, const GridFunction *gf,
77  const int *partitioning = NULL);
78 
79  /** @brief Construct a ParGridFunction on a given ParMesh, @a pmesh, reading
80  from an std::istream.
81 
82  In the process, a ParFiniteElementSpace and a FiniteElementCollection are
83  constructed. The new ParGridFunction assumes ownership of both. */
84  ParGridFunction(ParMesh *pmesh, std::istream &input);
85 
86  /// Copy assignment. Only the data of the base class Vector is copied.
87  /** It is assumed that this object and @a rhs use ParFiniteElementSpace%s
88  that have the same size.
89 
90  @note Defining this method overwrites the implicitly defined copy
91  assignemnt operator. */
93  { return operator=((const Vector &)rhs); }
94 
95  /// Assign constant values to the ParGridFunction data.
96  ParGridFunction &operator=(double value)
97  { GridFunction::operator=(value); return *this; }
98 
99  /// Copy the data from a Vector to the ParGridFunction data.
101  { GridFunction::operator=(v); return *this; }
102 
103  ParFiniteElementSpace *ParFESpace() const { return pfes; }
104 
105  virtual void Update();
106 
107  /// Associate a new FiniteElementSpace with the ParGridFunction.
108  /** The ParGridFunction is resized using the SetSize() method. The new space
109  @a f is expected to be a ParFiniteElementSpace. */
110  virtual void SetSpace(FiniteElementSpace *f);
111 
112  /// Associate a new parallel space with the ParGridFunction.
114 
115  /** @brief Make the ParGridFunction reference external data on a new
116  FiniteElementSpace. */
117  /** This method changes the FiniteElementSpace associated with the
118  ParGridFunction and sets the pointer @a v as external data in the
119  ParGridFunction. The new space @a f is expected to be a
120  ParFiniteElementSpace. */
121  virtual void MakeRef(FiniteElementSpace *f, double *v);
122 
123  /** @brief Make the ParGridFunction reference external data on a new
124  ParFiniteElementSpace. */
125  /** This method changes the ParFiniteElementSpace associated with the
126  ParGridFunction and sets the pointer @a v as external data in the
127  ParGridFunction. */
128  void MakeRef(ParFiniteElementSpace *f, double *v);
129 
130  /** @brief Make the ParGridFunction reference external data on a new
131  FiniteElementSpace. */
132  /** This method changes the FiniteElementSpace associated with the
133  ParGridFunction and sets the data of the Vector @a v (plus the @a
134  v_offset) as external data in the ParGridFunction. The new space @a f is
135  expected to be a ParFiniteElementSpace.
136  @note This version of the method will also perform bounds checks when
137  the build option MFEM_DEBUG is enabled. */
138  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
139 
140  /** @brief Make the ParGridFunction reference external data on a new
141  ParFiniteElementSpace. */
142  /** This method changes the ParFiniteElementSpace associated with the
143  ParGridFunction and sets the data of the Vector @a v (plus the
144  @a v_offset) as external data in the ParGridFunction.
145  @note This version of the method will also perform bounds checks when
146  the build option MFEM_DEBUG is enabled. */
147  void MakeRef(ParFiniteElementSpace *f, Vector &v, int v_offset);
148 
149  /** Set the grid function on (all) dofs from a given vector on the
150  true dofs, i.e. P tv. */
151  void Distribute(const Vector *tv);
152  void Distribute(const Vector &tv) { Distribute(&tv); }
153  void AddDistribute(double a, const Vector *tv);
154  void AddDistribute(double a, const Vector &tv) { AddDistribute(a, &tv); }
155 
156  /// Set the GridFunction from the given true-dof vector.
157  virtual void SetFromTrueDofs(const Vector &tv) { Distribute(tv); }
158 
159  /// Short semantic for Distribute()
161  { Distribute(&tv); return (*this); }
162 
164 
165  /// Returns the true dofs in a new HypreParVector
166  HypreParVector *GetTrueDofs() const;
167 
168  /// Returns the vector averaged on the true dofs.
169  void ParallelAverage(Vector &tv) const;
170 
171  /// Returns the vector averaged on the true dofs.
172  void ParallelAverage(HypreParVector &tv) const;
173 
174  /// Returns a new vector averaged on the true dofs.
176 
177  /// Returns the vector restricted to the true dofs.
178  void ParallelProject(Vector &tv) const;
179 
180  /// Returns the vector restricted to the true dofs.
181  void ParallelProject(HypreParVector &tv) const;
182 
183  /// Returns a new vector restricted to the true dofs.
185 
186  /// Returns the vector assembled on the true dofs.
187  void ParallelAssemble(Vector &tv) const;
188 
189  /// Returns the vector assembled on the true dofs.
190  void ParallelAssemble(HypreParVector &tv) const;
191 
192  /// Returns a new vector assembled on the true dofs.
194 
195  void ExchangeFaceNbrData();
197  const Vector &FaceNbrData() const { return face_nbr_data; }
198 
199  // Redefine to handle the case when i is a face-neighbor element
200  virtual double GetValue(int i, const IntegrationPoint &ip,
201  int vdim = 1) const;
203  { return GetValue(T.ElementNo, T.GetIntPoint()); }
204 
206  virtual void ProjectCoefficient(Coefficient &coeff);
207 
209  /** @brief Project a discontinuous vector coefficient as a grid function on
210  a continuous finite element space. The values in shared dofs are
211  determined from the element with maximal attribute. */
212  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
213 
214  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
215 
216  virtual void ProjectDiscCoefficient(VectorCoefficient &vcoeff, AvgType type);
217 
219 
220  // Only the values in the master are guaranteed to be correct!
222  Array<int> &attr)
223  { ProjectBdrCoefficient(NULL, &vcoeff, attr); }
224 
225  // Only the values in the master are guaranteed to be correct!
226  virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr)
227  { ProjectBdrCoefficient(coeff, NULL, attr); }
228 
229  // Only the values in the master are guaranteed to be correct!
230  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
231  Array<int> &bdr_attr);
232 
233  virtual double ComputeL1Error(Coefficient *exsol[],
234  const IntegrationRule *irs[] = NULL) const
235  {
237  *exsol, NULL, 1, NULL, irs), pfes->GetComm());
238  }
239 
240  virtual double ComputeL1Error(Coefficient &exsol,
241  const IntegrationRule *irs[] = NULL) const
242  { return ComputeLpError(1.0, exsol, NULL, irs); }
243 
244  virtual double ComputeL1Error(VectorCoefficient &exsol,
245  const IntegrationRule *irs[] = NULL) const
246  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
247 
248  virtual double ComputeL2Error(Coefficient *exsol[],
249  const IntegrationRule *irs[] = NULL) const
250  {
251  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs),
252  pfes->GetComm());
253  }
254 
255  virtual double ComputeL2Error(Coefficient &exsol,
256  const IntegrationRule *irs[] = NULL) const
257  { return ComputeLpError(2.0, exsol, NULL, irs); }
258 
259  virtual double ComputeL2Error(VectorCoefficient &exsol,
260  const IntegrationRule *irs[] = NULL,
261  Array<int> *elems = NULL) const
262  {
263  return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
264  pfes->GetComm());
265  }
266 
267  virtual double ComputeMaxError(Coefficient *exsol[],
268  const IntegrationRule *irs[] = NULL) const
269  {
270  return GlobalLpNorm(infinity(),
271  GridFunction::ComputeMaxError(exsol, irs),
272  pfes->GetComm());
273  }
274 
275  virtual double ComputeMaxError(Coefficient &exsol,
276  const IntegrationRule *irs[] = NULL) const
277  {
278  return ComputeLpError(infinity(), exsol, NULL, irs);
279  }
280 
281  virtual double ComputeMaxError(VectorCoefficient &exsol,
282  const IntegrationRule *irs[] = NULL) const
283  {
284  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
285  }
286 
287  virtual double ComputeLpError(const double p, Coefficient &exsol,
288  Coefficient *weight = NULL,
289  const IntegrationRule *irs[] = NULL) const
290  {
292  p, exsol, weight, irs), pfes->GetComm());
293  }
294 
295  /** When given a vector weight, compute the pointwise (scalar) error as the
296  dot product of the vector error with the vector weight. Otherwise, the
297  scalar error is the l_2 norm of the vector error. */
298  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
299  Coefficient *weight = NULL,
300  VectorCoefficient *v_weight = NULL,
301  const IntegrationRule *irs[] = NULL) const
302  {
304  p, exsol, weight, v_weight, irs), pfes->GetComm());
305  }
306 
307  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
308  GridFunction &flux,
309  int wcoef = 1, int subdomain = -1);
310 
311  /** Save the local portion of the ParGridFunction. It differs from the
312  serial GridFunction::Save in that it takes into account the signs of
313  the local dofs. */
314  virtual void Save(std::ostream &out) const;
315 
316  /// Merge the local grid functions
317  void SaveAsOne(std::ostream &out = mfem::out);
318 
319  virtual ~ParGridFunction() { }
320 };
321 
322 
323 /** Performs a global L2 projection (through a HypreBoomerAMG solve) of flux
324  from supplied discontinuous space into supplied smooth (continuous, or at
325  least conforming) space, and computes the Lp norms of the differences
326  between them on each element. This is one approach to handling conforming
327  and non-conforming elements in parallel. Returns the total error estimate. */
328 double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator,
329  const ParGridFunction &x,
330  ParFiniteElementSpace &smooth_flux_fes,
331  ParFiniteElementSpace &flux_fes,
332  Vector &errors, int norm_p = 2, double solver_tol = 1e-12,
333  int solver_max_it = 200);
334 
335 }
336 
337 #endif // MFEM_USE_MPI
338 
339 #endif
virtual double ComputeLpError(const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:298
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:85
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Memory< double > data
Definition: vector.hpp:52
void AddDistribute(double a, const Vector &tv)
Definition: pgridfunc.hpp:154
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
Definition: pgridfunc.hpp:48
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
Definition: pgridfunc.cpp:196
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:306
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
Definition: pgridfunc.cpp:179
virtual ~ParGridFunction()
Definition: pgridfunc.hpp:319
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1713
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:160
void SaveAsOne(std::ostream &out=mfem::out)
Merge the local grid functions.
Definition: pgridfunc.cpp:504
ParGridFunction(ParFiniteElementSpace *pf)
Definition: pgridfunc.hpp:51
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: pgridfunc.cpp:254
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:488
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Definition: pgridfunc.cpp:313
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:86
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:292
const IntegrationPoint & GetIntPoint()
Definition: eltrans.hpp:55
void Distribute(const Vector &tv)
Definition: pgridfunc.hpp:152
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:221
double GetValue(ElementTransformation &T)
Definition: pgridfunc.hpp:202
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:101
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
ParGridFunction & operator=(double value)
Assign constant values to the ParGridFunction data.
Definition: pgridfunc.hpp:96
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:157
const Vector & FaceNbrData() const
Definition: pgridfunc.hpp:197
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:80
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:244
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction. If all dofs are true, then tv will be set to point to th...
Definition: gridfunc.cpp:323
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:287
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: gridfunc.hpp:106
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2305
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
Definition: pgridfunc.hpp:100
Vector & FaceNbrData()
Definition: pgridfunc.hpp:196
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:233
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:636
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:23
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: pgridfunc.hpp:92
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:255
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:400
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
virtual double ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:267
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:321
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:709
void Distribute(const Vector *tv)
Definition: pgridfunc.cpp:131
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:275
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:281
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:226
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:248
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:444
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
Definition: pgridfunc.cpp:162
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:141
void AddDistribute(double a, const Vector *tv)
Definition: pgridfunc.cpp:136
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1582
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:240
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2199
virtual double ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const
Definition: pgridfunc.hpp:259
Vector data type.
Definition: vector.hpp:48
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:64
Class for parallel meshes.
Definition: pmesh.hpp:32
ParGridFunction(ParFiniteElementSpace *pf, double *data)
Construct a ParGridFunction using previously allocated array data.
Definition: pgridfunc.hpp:60
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: pgridfunc.cpp:672
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:275
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:103