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