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