MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
gridfunc.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_GRIDFUNC
13 #define MFEM_GRIDFUNC
14 
15 #include "../config/config.hpp"
16 #include "fespace.hpp"
17 #include "coefficient.hpp"
18 #include "bilininteg.hpp"
19 #include <limits>
20 #include <ostream>
21 #include <string>
22 
23 namespace mfem
24 {
25 
26 /// Class for grid function - Vector with associated FE space.
27 class GridFunction : public Vector
28 {
29 protected:
30  /// FE space on which grid function lives.
32 
33  /// Used when the grid function is read from a file
35 
36  long sequence; // see FiniteElementSpace::sequence, Mesh::sequence
37 
38  void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[]);
39 
41 
42  // Project the delta coefficient without scaling and return the (local)
43  // integral of the projection.
44  void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff,
45  double &integral);
46 
47  // Sum fluxes to vertices and count element contributions
49  GridFunction &flux,
50  Array<int>& counts,
51  int wcoef,
52  int subdomain);
53 
54  /** Project a discontinuous vector coefficient in a continuous space and
55  return in dof_attr the maximal attribute of the elements containing each
56  degree of freedom. */
57  void ProjectDiscCoefficient(VectorCoefficient &coeff, Array<int> &dof_attr);
58 
59  void Destroy();
60 
61 public:
62 
63  GridFunction() { fes = NULL; fec = NULL; sequence = 0; }
64 
65  /// Copy constructor.
67  : Vector(orig), fes(orig.fes), fec(NULL), sequence(orig.sequence) { }
68 
69  /// Construct a GridFunction associated with the FiniteElementSpace @a *f.
70  GridFunction(FiniteElementSpace *f) : Vector(f->GetVSize())
71  { fes = f; fec = NULL; sequence = f->GetSequence(); }
72 
73  /// Construct a GridFunction using previously allocated array @a data.
74  /** The GridFunction does not assume ownership of @a data which is assumed to
75  be of size at least `f->GetVSize()`. Similar to the Vector constructor
76  for externally allocated array, the pointer @a data can be NULL. The data
77  array can be replaced later using the method SetData().
78  */
79  GridFunction(FiniteElementSpace *f, double *data) : Vector(data, f->GetVSize())
80  { fes = f; fec = NULL; sequence = f->GetSequence(); }
81 
82  /// Construct a GridFunction on the given Mesh, using the data from @a input.
83  /** The content of @a input should be in the format created by the method
84  Save(). The reconstructed FiniteElementSpace and FiniteElementCollection
85  are owned by the GridFunction. */
86  GridFunction(Mesh *m, std::istream &input);
87 
88  GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces);
89 
90  /// Make the GridFunction the owner of 'fec' and 'fes'
91  void MakeOwner(FiniteElementCollection *_fec) { fec = _fec; }
92 
94 
95  int VectorDim() const;
96 
97  /// @brief Extract the true-dofs from the GridFunction. If all dofs are true,
98  /// then `tv` will be set to point to the data of `*this`.
99  void GetTrueDofs(Vector &tv) const;
100 
101  /// Set the GridFunction from the given true-dof vector.
102  virtual void SetFromTrueDofs(const Vector &tv);
103 
104  /// Returns the values in the vertices of i'th element for dimension vdim.
105  void GetNodalValues(int i, Array<double> &nval, int vdim = 1) const;
106 
107  virtual double GetValue(int i, const IntegrationPoint &ip,
108  int vdim = 1) const;
109 
110  void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const;
111 
112  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
113  int vdim = 1) const;
114 
115  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
116  DenseMatrix &tr, int vdim = 1) const;
117 
118  int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals,
119  DenseMatrix &tr, int vdim = 1) const;
120 
122  DenseMatrix &vals) const;
123 
124  void GetVectorValues(int i, const IntegrationRule &ir,
125  DenseMatrix &vals, DenseMatrix &tr) const;
126 
127  int GetFaceVectorValues(int i, int side, const IntegrationRule &ir,
128  DenseMatrix &vals, DenseMatrix &tr) const;
129 
130  void GetValuesFrom(GridFunction &);
131 
133 
134  void GetVectorFieldValues(int i, const IntegrationRule &ir,
135  DenseMatrix &vals,
136  DenseMatrix &tr, int comp = 0) const;
137 
138  /// For a vector grid function, makes sure that the ordering is byNODES.
139  void ReorderByNodes();
140 
141  /// Return the values as a vector on mesh vertices for dimension vdim.
142  void GetNodalValues(Vector &nval, int vdim = 1) const;
143 
144  void GetVectorFieldNodalValues(Vector &val, int comp) const;
145 
146  void ProjectVectorFieldOn(GridFunction &vec_field, int comp = 0);
147 
148  void GetDerivative(int comp, int der_comp, GridFunction &der);
149 
151 
152  void GetCurl(ElementTransformation &tr, Vector &curl);
153 
154  void GetGradient(ElementTransformation &tr, Vector &grad);
155 
156  void GetGradients(const int elem, const IntegrationRule &ir,
157  DenseMatrix &grad);
158 
160 
161  /** Compute \f$ (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) \f$,
162  where \f$ \psi_i \f$ are the basis functions for the FE space of avgs.
163  Both FE spaces should be scalar and on the same mesh. */
164  void GetElementAverages(GridFunction &avgs);
165 
166  /** Impose the given bounds on the function's DOFs while preserving its local
167  * integral (described in terms of the given weights) on the i'th element
168  * through SLBPQ optimization.
169  * Intended to be used for discontinuous FE functions. */
170  void ImposeBounds(int i, const Vector &weights,
171  const Vector &_lo, const Vector &_hi);
172  void ImposeBounds(int i, const Vector &weights,
173  double _min = 0.0, double _max = std::numeric_limits<double>::infinity());
174 
175  /** Project the given 'src' GridFunction to 'this' GridFunction, both of
176  which must be on the same mesh. The current implementation assumes that
177  all element use the same projection matrix. */
178  void ProjectGridFunction(const GridFunction &src);
179 
180  virtual void ProjectCoefficient(Coefficient &coeff);
181 
182  // call fes -> BuildDofToArrays() before using this projection
183  void ProjectCoefficient(Coefficient &coeff, Array<int> &dofs, int vd = 0);
184 
186 
187  // call fes -> BuildDofToArrays() before using this projection
188  void ProjectCoefficient(VectorCoefficient &vcoeff, Array<int> &dofs);
189 
190  void ProjectCoefficient(Coefficient *coeff[]);
191 
192  /** @brief Project a discontinuous vector coefficient as a grid function on
193  a continuous finite element space. The values in shared dofs are
194  determined from the element with maximal attribute. */
195  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
196 
198  /** @brief Projects a discontinuous coefficient so that the values in shared
199  vdofs are computed by taking an average of the possible values. */
200  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
201 
202 protected:
203  /** @brief Accumulates (depending on @a type) the values of @a coeff at all
204  shared vdofs and counts in how many zones each vdof appears. */
205  void AccumulateAndCountZones(Coefficient &coeff, AvgType type,
206  Array<int> &zones_per_vdof);
207 
208  // Complete the computation of averages; called e.g. after
209  // AccumulateAndCountZones().
210  void ComputeMeans(AvgType type, Array<int> &zones_per_vdof);
211 
212 public:
214  {
215  Coefficient *coeff_p = &coeff;
216  ProjectBdrCoefficient(&coeff_p, attr);
217  }
218 
219  void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
220 
221  /** Project the normal component of the given VectorCoefficient on
222  the boundary. Only boundary attributes that are marked in
223  'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
225  Array<int> &bdr_attr);
226 
227  /** Project the tangential components of the given VectorCoefficient on
228  the boundary. Only boundary attributes that are marked in
229  'bdr_attr' are projected. Assumes ND-type VectorFE GridFunction. */
231  Array<int> &bdr_attr);
232 
233  virtual double ComputeL2Error(Coefficient &exsol,
234  const IntegrationRule *irs[] = NULL) const
235  { return ComputeLpError(2.0, exsol, NULL, irs); }
236 
237  virtual double ComputeL2Error(Coefficient *exsol[],
238  const IntegrationRule *irs[] = NULL) const;
239 
240  virtual double ComputeL2Error(VectorCoefficient &exsol,
241  const IntegrationRule *irs[] = NULL,
242  Array<int> *elems = NULL) const;
243 
244  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
245  Coefficient *ell_coef, double Nu,
246  int norm_type) const;
247 
248  virtual double ComputeMaxError(Coefficient &exsol,
249  const IntegrationRule *irs[] = NULL) const
250  {
251  return ComputeLpError(std::numeric_limits<double>::infinity(),
252  exsol, NULL, irs);
253  }
254 
255  virtual double ComputeMaxError(Coefficient *exsol[],
256  const IntegrationRule *irs[] = NULL) const;
257 
258  virtual double ComputeMaxError(VectorCoefficient &exsol,
259  const IntegrationRule *irs[] = NULL) const
260  {
261  return ComputeLpError(std::numeric_limits<double>::infinity(),
262  exsol, NULL, NULL, irs);
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 ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
270  int norm_type, Array<int> *elems = NULL,
271  const IntegrationRule *irs[] = NULL) const;
272 
273  virtual double ComputeL1Error(VectorCoefficient &exsol,
274  const IntegrationRule *irs[] = NULL) const
275  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
276 
277  virtual double ComputeLpError(const double p, Coefficient &exsol,
278  Coefficient *weight = NULL,
279  const IntegrationRule *irs[] = NULL) const;
280 
281  /** When given a vector weight, compute the pointwise (scalar) error as the
282  dot product of the vector error with the vector weight. Otherwise, the
283  scalar error is the l_2 norm of the vector error. */
284  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
285  Coefficient *weight = NULL,
286  VectorCoefficient *v_weight = NULL,
287  const IntegrationRule *irs[] = NULL) const;
288 
289  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
290  GridFunction &flux,
291  int wcoef = 1, int subdomain = -1);
292 
293  /// Redefine '=' for GridFunction = constant.
294  GridFunction &operator=(double value);
295 
296  /// Copy the data from @a v.
297  /** The size of @a v must be equal to the size of the FiniteElementSpace
298  @a fes. */
299  GridFunction &operator=(const Vector &v);
300 
301  /// Copy the data from @a v.
302  /** The GridFunctions @a v and @a *this must have FiniteElementSpaces with
303  the same size. */
305 
306  /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
307  virtual void Update();
308 
310  const FiniteElementSpace *FESpace() const { return fes; }
311 
312  /// Associate a new FiniteElementSpace with the GridFunction.
313  /** The GridFunction is resized using the SetSize() method. */
314  virtual void SetSpace(FiniteElementSpace *f);
315 
316  /** @brief Make the GridFunction reference external data on a new
317  FiniteElementSpace. */
318  /** This method changes the FiniteElementSpace associated with the
319  GridFunction and sets the pointer @a v as external data in the
320  GridFunction. */
321  virtual void MakeRef(FiniteElementSpace *f, double *v);
322 
323  /** @brief Make the GridFunction reference external data on a new
324  FiniteElementSpace. */
325  /** This method changes the FiniteElementSpace associated with the
326  GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
327  as external data in the GridFunction.
328  @note This version of the method will also perform bounds checks when
329  the build option MFEM_DEBUG is enabled. */
330  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
331 
332  /// Save the GridFunction to an output stream.
333  virtual void Save(std::ostream &out) const;
334 
335  /** Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be
336  called first. The parameter ref > 0 must match the one used in
337  Mesh::PrintVTK. */
338  void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
339 
340  void SaveSTL(std::ostream &out, int TimesToRefine = 1);
341 
342  /// Destroys grid function.
343  virtual ~GridFunction() { Destroy(); }
344 };
345 
346 
347 /** Overload operator<< for std::ostream and GridFunction; valid also for the
348  derived class ParGridFunction */
349 std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
350 
351 
352 /** @brief Class representing a function through its values (scalar or vector)
353  at quadrature points. */
355 {
356 protected:
357  QuadratureSpace *qspace; ///< Associated QuadratureSpace
358  int vdim; ///< Vector dimension
359  bool own_qspace; ///< QuadratureSpace ownership flag
360 
361 public:
362  /// Create an empty QuadratureFunction.
363  /** The object can be initialized later using the SetSpace() methods. */
365  : qspace(NULL), vdim(0), own_qspace(false) { }
366 
367  /// Create a QuadratureFunction based on the given QuadratureSpace.
368  /** The QuadratureFunction does not assume ownership of the QuadratureSpace.
369  @note The Vector data is not initialized. */
370  QuadratureFunction(QuadratureSpace *qspace_, int vdim_ = 1)
371  : Vector(vdim_*qspace_->GetSize()),
372  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
373 
374  /** @brief Create a QuadratureFunction based on the given QuadratureSpace,
375  using the external data, @a qf_data. */
376  /** The QuadratureFunction does not assume ownership of neither the
377  QuadratureSpace nor the external data. */
378  QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_ = 1)
379  : Vector(qf_data, vdim_*qspace_->GetSize()),
380  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
381 
382  /// Read a QuadratureFunction from the stream @a in.
383  /** The QuadratureFunction assumes ownership of the read QuadratureSpace. */
384  QuadratureFunction(Mesh *mesh, std::istream &in);
385 
386  virtual ~QuadratureFunction() { if (own_qspace) { delete qspace; } }
387 
388  /// Get the associated QuadratureSpace.
389  QuadratureSpace *GetSpace() const { return qspace; }
390 
391  /// Change the QuadratureSpace and optionally the vector dimension.
392  /** If the new QuadratureSpace is different from the current one, the
393  QuadratureFunction will not assume ownership of the new space; otherwise,
394  the ownership flag remains the same.
395 
396  If the new vector dimension @a vdim_ < 0, the vector dimension remains
397  the same.
398 
399  The data size is updated by calling Vector::SetSize(). */
400  inline void SetSpace(QuadratureSpace *qspace_, int vdim_ = -1);
401 
402  /** @brief Change the QuadratureSpace, the data array, and optionally the
403  vector dimension. */
404  /** If the new QuadratureSpace is different from the current one, the
405  QuadratureFunction will not assume ownership of the new space; otherwise,
406  the ownership flag remains the same.
407 
408  If the new vector dimension @a vdim_ < 0, the vector dimension remains
409  the same.
410 
411  The data array is replaced by calling Vector::NewDataAndSize(). */
412  inline void SetSpace(QuadratureSpace *qspace_, double *qf_data,
413  int vdim_ = -1);
414 
415  /// Get the vector dimension.
416  int GetVDim() const { return vdim; }
417 
418  /// Set the vector dimension, updating the size by calling Vector::SetSize().
419  void SetVDim(int vdim_)
420  { vdim = vdim_; SetSize(vdim*qspace->GetSize()); }
421 
422  /// Get the QuadratureSpace ownership flag.
423  bool OwnsSpace() { return own_qspace; }
424 
425  /// Set the QuadratureSpace ownership flag.
426  void SetOwnsSpace(bool own) { own_qspace = own; }
427 
428  /// Get the IntegrationRule associated with mesh element @a idx.
430  { return qspace->GetElementIntRule(idx); }
431 
432  /// Return all values associated with mesh element @a idx in a Vector.
433  /** The result is stored in the Vector @a values as a reference to the
434  global values.
435 
436  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
437  `i`-th vector component at the `j`-th quadrature point.
438  */
439  inline void GetElementValues(int idx, Vector &values);
440 
441  /// Return all values associated with mesh element @a idx in a DenseMatrix.
442  /** The result is stored in the DenseMatrix @a values as a reference to the
443  global values.
444 
445  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
446  `i`-th vector component at the `j`-th quadrature point.
447  */
448  inline void GetElementValues(int idx, DenseMatrix &values);
449 
450  /// Write the QuadratureFunction to the stream @a out.
451  void Save(std::ostream &out) const;
452 };
453 
454 /// Overload operator<< for std::ostream and QuadratureFunction.
455 std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
456 
457 
458 double ZZErrorEstimator(BilinearFormIntegrator &blfi,
459  GridFunction &u,
460  GridFunction &flux,
461  Vector &error_estimates,
462  Array<int> *aniso_flags = NULL,
463  int with_subdomains = 1);
464 
465 /// Compute the Lp distance between two grid functions on the given element.
466 double ComputeElementLpDistance(double p, int i,
467  GridFunction& gf1, GridFunction& gf2);
468 
469 
470 /// Class used for extruding scalar GridFunctions
472 {
473 private:
474  int n;
475  Mesh *mesh_in;
476  Coefficient &sol_in;
477 public:
479  : n(_n), mesh_in(m), sol_in(s) { }
480  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
481  virtual ~ExtrudeCoefficient() { }
482 };
483 
484 /// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
485 GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
486  GridFunction *sol, const int ny);
487 
488 
489 // Inline methods
490 
491 inline void QuadratureFunction::SetSpace(QuadratureSpace *qspace_, int vdim_)
492 {
493  if (qspace_ != qspace)
494  {
495  if (own_qspace) { delete qspace; }
496  qspace = qspace_;
497  own_qspace = false;
498  }
499  vdim = (vdim_ < 0) ? vdim : vdim_;
501 }
502 
504  double *qf_data, int vdim_)
505 {
506  if (qspace_ != qspace)
507  {
508  if (own_qspace) { delete qspace; }
509  qspace = qspace_;
510  own_qspace = false;
511  }
512  vdim = (vdim_ < 0) ? vdim : vdim_;
513  NewDataAndSize(qf_data, vdim*qspace->GetSize());
514 }
515 
516 inline void QuadratureFunction::GetElementValues(int idx, Vector &values)
517 {
518  const int s_offset = qspace->element_offsets[idx];
519  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
520  values.NewDataAndSize(data + vdim*s_offset, vdim*sl_size);
521 }
522 
524 {
525  const int s_offset = qspace->element_offsets[idx];
526  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
527  values.Reset(data + vdim*s_offset, vdim, sl_size);
528 }
529 
530 } // namespace mfem
531 
532 #endif
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:2351
int GetVDim() const
Get the vector dimension.
Definition: gridfunc.hpp:416
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.hpp:343
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8533
QuadratureFunction(QuadratureSpace *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace.
Definition: gridfunc.hpp:370
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:101
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:471
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:233
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1232
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains)
Definition: gridfunc.cpp:2502
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:370
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:320
GridFunction * Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d, GridFunction *sol, const int ny)
Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
Definition: gridfunc.cpp:2656
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of &#39;fec&#39; and &#39;fes&#39;.
Definition: gridfunc.hpp:91
int GetSize()
Return the total number of quadrature points.
Definition: fespace.hpp:415
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh)
Definition: gridfunc.cpp:858
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: gridfunc.cpp:246
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1424
bool own_qspace
QuadratureSpace ownership flag.
Definition: gridfunc.hpp:359
void GetBdrValuesFrom(GridFunction &)
Definition: gridfunc.cpp:621
void GetGradient(ElementTransformation &tr, Vector &grad)
Definition: gridfunc.cpp:966
int VectorDim() const
Definition: gridfunc.cpp:261
virtual ~QuadratureFunction()
Definition: gridfunc.hpp:386
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient.
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
Definition: gridfunc.cpp:1087
int vdim
Vector dimension.
Definition: gridfunc.hpp:358
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:258
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
Definition: gridfunc.hpp:491
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:1209
bool OwnsSpace()
Get the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:423
void GetNodalValues(int i, Array< double > &nval, int vdim=1) const
Returns the values in the vertices of i&#39;th element for dimension vdim.
Definition: gridfunc.cpp:314
QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace, using the external data...
Definition: gridfunc.hpp:378
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:727
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2259
const FiniteElementSpace * FESpace() const
Definition: gridfunc.hpp:310
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: gridfunc.cpp:2646
GridFunction(FiniteElementSpace *f, double *data)
Construct a GridFunction using previously allocated array data.
Definition: gridfunc.hpp:79
void GetDerivative(int comp, int der_comp, GridFunction &der)
Definition: gridfunc.cpp:798
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1549
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:199
double GetDivergence(ElementTransformation &tr)
Definition: gridfunc.cpp:877
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:265
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:71
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:658
void SetOwnsSpace(bool own)
Set the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:426
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:1767
GridFunction(const GridFunction &orig)
Copy constructor.
Definition: gridfunc.hpp:66
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:389
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:273
QuadratureFunction()
Create an empty QuadratureFunction.
Definition: gridfunc.hpp:364
void AccumulateAndCountZones(Coefficient &coeff, AvgType type, Array< int > &zones_per_vdof)
Accumulates (depending on type) the values of coeff at all shared vdofs and counts in how many zones ...
Definition: gridfunc.cpp:1171
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:181
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition: gridfunc.hpp:419
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
void GetCurl(ElementTransformation &tr, Vector &curl)
Definition: gridfunc.cpp:913
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:404
FiniteElementCollection * OwnFEC()
Definition: gridfunc.hpp:93
virtual ~ExtrudeCoefficient()
Definition: gridfunc.hpp:481
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:309
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:418
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2078
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Definition: gridfunc.hpp:429
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
FiniteElementSpace * fes
FE space on which grid function lives.
Definition: gridfunc.hpp:31
void ProjectGridFunction(const GridFunction &src)
Definition: gridfunc.cpp:1054
QuadratureSpace * qspace
Associated QuadratureSpace.
Definition: gridfunc.hpp:357
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:248
FiniteElementCollection * fec
Used when the grid function is read from a file.
Definition: gridfunc.hpp:34
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:146
double ComputeElementLpDistance(double p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
Definition: gridfunc.cpp:2589
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:541
Class for integration point with weight.
Definition: intrules.hpp:25
void SaveSTL(std::ostream &out, int TimesToRefine=1)
Definition: gridfunc.cpp:2371
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:437
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: gridfunc.hpp:516
void GetElementAverages(GridFunction &avgs)
Definition: gridfunc.cpp:1024
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1620
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:487
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1295
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:173
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:700
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:1972
ExtrudeCoefficient(Mesh *m, Coefficient &s, int _n)
Definition: gridfunc.hpp:478
GridFunction & operator=(double value)
Redefine &#39;=&#39; for GridFunction = constant.
Definition: gridfunc.cpp:2234
Vector data type.
Definition: vector.hpp:41
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad)
Definition: gridfunc.cpp:985
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Definition: gridfunc.hpp:70
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:386
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
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:376
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Definition: gridfunc.cpp:1010
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: gridfunc.cpp:297
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:354
void GetValuesFrom(GridFunction &)
Definition: gridfunc.cpp:584
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
Definition: gridfunc.cpp:2486
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
Definition: gridfunc.cpp:2274
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Definition: gridfunc.hpp:213
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:355
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:758