MFEM  v3.3
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  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  /** Project a discontinuous vector coefficient as a grid function on a
193  continuous finite element space. The values in shared dofs are determined
194  from the element with maximal attribute. */
196 
198  {
199  Coefficient *coeff_p = &coeff;
200  ProjectBdrCoefficient(&coeff_p, attr);
201  }
202 
203  void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
204 
205  /** Project the normal component of the given VectorCoefficient on
206  the boundary. Only boundary attributes that are marked in
207  'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
209  Array<int> &bdr_attr);
210 
211  /** Project the tangential components of the given VectorCoefficient on
212  the boundary. Only boundary attributes that are marked in
213  'bdr_attr' are projected. Assumes ND-type VectorFE GridFunction. */
215  Array<int> &bdr_attr);
216 
218  const IntegrationRule *irs[] = NULL) const
219  { return ComputeLpError(2.0, exsol, NULL, irs); }
220 
221  double ComputeL2Error(Coefficient *exsol[],
222  const IntegrationRule *irs[] = NULL) const;
223 
224  double ComputeL2Error(VectorCoefficient &exsol,
225  const IntegrationRule *irs[] = NULL,
226  Array<int> *elems = NULL) const;
227 
228  double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
229  Coefficient *ell_coef, double Nu,
230  int norm_type) const;
231 
233  const IntegrationRule *irs[] = NULL) const
234  {
235  return ComputeLpError(std::numeric_limits<double>::infinity(),
236  exsol, NULL, irs);
237  }
238 
239  double ComputeMaxError(Coefficient *exsol[],
240  const IntegrationRule *irs[] = NULL) const;
241 
243  const IntegrationRule *irs[] = NULL) const
244  {
245  return ComputeLpError(std::numeric_limits<double>::infinity(),
246  exsol, NULL, NULL, irs);
247  }
248 
250  const IntegrationRule *irs[] = NULL) const
251  { return ComputeLpError(1.0, exsol, NULL, irs); }
252 
253  double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
254  int norm_type, Array<int> *elems = NULL,
255  const IntegrationRule *irs[] = NULL) const;
256 
258  const IntegrationRule *irs[] = NULL) const
259  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
260 
261  double ComputeLpError(const double p, Coefficient &exsol,
262  Coefficient *weight = NULL,
263  const IntegrationRule *irs[] = NULL) const;
264 
265  /** When given a vector weight, compute the pointwise (scalar) error as the
266  dot product of the vector error with the vector weight. Otherwise, the
267  scalar error is the l_2 norm of the vector error. */
268  double ComputeLpError(const double p, VectorCoefficient &exsol,
269  Coefficient *weight = NULL,
270  VectorCoefficient *v_weight = NULL,
271  const IntegrationRule *irs[] = NULL) const;
272 
273  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
274  GridFunction &flux,
275  int wcoef = 1, int subdomain = -1);
276 
277  /// Redefine '=' for GridFunction = constant.
278  GridFunction &operator=(double value);
279 
280  /// Copy the data from @a v.
281  /** The size of @a v must be equal to the size of the FiniteElementSpace
282  @a fes. */
283  GridFunction &operator=(const Vector &v);
284 
285  /// Copy the data from @a v.
286  /** The GridFunctions @a v and @a *this must have FiniteElementSpaces with
287  the same size. */
289 
290  /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
291  void Update();
292 
294  const FiniteElementSpace *FESpace() const { return fes; }
295 
296  void SetSpace(FiniteElementSpace *f);
297 
298  /** @brief Make the GridFunction reference external data on a new
299  FiniteElementSpace. */
300  /** This method changes the FiniteElementSpace associated with the
301  GridFunction and sets the pointer @a v as external data in the
302  GridFunction. */
303  void MakeRef(FiniteElementSpace *f, double *v);
304 
305  /** @brief Make the GridFunction reference external data on a new
306  FiniteElementSpace. */
307  /** This method changes the FiniteElementSpace associated with the
308  GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
309  as external data in the GridFunction.
310  @note This version of the method will also perform bounds checks when
311  the build option MFEM_DEBUG is enabled. */
312  void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
313 
314  /// Save the GridFunction to an output stream.
315  virtual void Save(std::ostream &out) const;
316 
317  /** Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be
318  called first. The parameter ref > 0 must match the one used in
319  Mesh::PrintVTK. */
320  void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
321 
322  void SaveSTL(std::ostream &out, int TimesToRefine = 1);
323 
324  /// Destroys grid function.
325  virtual ~GridFunction() { Destroy(); }
326 };
327 
328 
329 /** Overload operator<< for std::ostream and GridFunction; valid also for the
330  derived class ParGridFunction */
331 std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
332 
333 
334 /** @brief Class representing a function through its values (scalar or vector)
335  at quadrature points. */
337 {
338 protected:
339  QuadratureSpace *qspace; ///< Associated QuadratureSpace
340  int vdim; ///< Vector dimension
341  bool own_qspace; ///< QuadratureSpace ownership flag
342 
343 public:
344  /// Create an empty QuadratureFunction.
345  /** The object can be initialized later using the SetSpace() methods. */
347  : qspace(NULL), vdim(0), own_qspace(false) { }
348 
349  /// Create a QuadratureFunction based on the given QuadratureSpace.
350  /** The QuadratureFunction does not assume ownership of the QuadratureSpace.
351  @note The Vector data is not initialized. */
352  QuadratureFunction(QuadratureSpace *qspace_, int vdim_ = 1)
353  : Vector(vdim_*qspace_->GetSize()),
354  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
355 
356  /** @brief Create a QuadratureFunction based on the given QuadratureSpace,
357  using the external data, @a qf_data. */
358  /** The QuadratureFunction does not assume ownership of neither the
359  QuadratureSpace nor the external data. */
360  QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_ = 1)
361  : Vector(qf_data, vdim_*qspace_->GetSize()),
362  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
363 
364  /// Read a QuadratureFunction from the stream @a in.
365  /** The QuadratureFunction assumes ownership of the read QuadratureSpace. */
366  QuadratureFunction(Mesh *mesh, std::istream &in);
367 
368  virtual ~QuadratureFunction() { if (own_qspace) { delete qspace; } }
369 
370  /// Get the associated QuadratureSpace.
371  QuadratureSpace *GetSpace() const { return qspace; }
372 
373  /// Change the QuadratureSpace and optionally the vector dimension.
374  /** If the new QuadratureSpace is different from the current one, the
375  QuadratureFunction will not assume ownership of the new space; otherwise,
376  the ownership flag remains the same.
377 
378  If the new vector dimension @a vdim_ < 0, the vector dimension remains
379  the same.
380 
381  The data size is updated by calling Vector::SetSize(). */
382  inline void SetSpace(QuadratureSpace *qspace_, int vdim_ = -1);
383 
384  /** @brief Change the QuadratureSpace, the data array, and optionally the
385  vector dimension. */
386  /** If the new QuadratureSpace is different from the current one, the
387  QuadratureFunction will not assume ownership of the new space; otherwise,
388  the ownership flag remains the same.
389 
390  If the new vector dimension @a vdim_ < 0, the vector dimension remains
391  the same.
392 
393  The data array is replaced by calling Vector::NewDataAndSize(). */
394  inline void SetSpace(QuadratureSpace *qspace_, double *qf_data,
395  int vdim_ = -1);
396 
397  /// Get the vector dimension.
398  int GetVDim() const { return vdim; }
399 
400  /// Set the vector dimension, updating the size by calling Vector::SetSize().
401  void SetVDim(int vdim_)
402  { vdim = vdim_; SetSize(vdim*qspace->GetSize()); }
403 
404  /// Get the QuadratureSpace ownership flag.
405  bool OwnsSpace() { return own_qspace; }
406 
407  /// Set the QuadratureSpace ownership flag.
408  void SetOwnsSpace(bool own) { own_qspace = own; }
409 
410  /// Get the IntegrationRule associated with mesh element @a idx.
412  { return qspace->GetElementIntRule(idx); }
413 
414  /// Return all values associated with mesh element @a idx in a Vector.
415  /** The result is stored in the Vector @a values as a reference to the
416  global values.
417 
418  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
419  `i`-th vector component at the `j`-th quadrature point.
420  */
421  inline void GetElementValues(int idx, Vector &values);
422 
423  /// Return all values associated with mesh element @a idx in a DenseMatrix.
424  /** The result is stored in the DenseMatrix @a values as a reference to the
425  global values.
426 
427  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
428  `i`-th vector component at the `j`-th quadrature point.
429  */
430  inline void GetElementValues(int idx, DenseMatrix &values);
431 
432  /// Write the QuadratureFunction to the stream @a out.
433  void Save(std::ostream &out) const;
434 };
435 
436 /// Overload operator<< for std::ostream and QuadratureFunction.
437 std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
438 
439 
440 double ZZErrorEstimator(BilinearFormIntegrator &blfi,
441  GridFunction &u,
442  GridFunction &flux,
443  Vector &error_estimates,
444  Array<int> *aniso_flags = NULL,
445  int with_subdomains = 1);
446 
447 /// Compute the Lp distance between two grid functions on the given element.
448 double ComputeElementLpDistance(double p, int i,
449  GridFunction& gf1, GridFunction& gf2);
450 
451 
452 /// Class used for extruding scalar GridFunctions
454 {
455 private:
456  int n;
457  Mesh *mesh_in;
458  Coefficient &sol_in;
459 public:
461  : n(_n), mesh_in(m), sol_in(s) { }
462  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
463  virtual ~ExtrudeCoefficient() { }
464 };
465 
466 /// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
467 GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
468  GridFunction *sol, const int ny);
469 
470 
471 // Inline methods
472 
473 inline void QuadratureFunction::SetSpace(QuadratureSpace *qspace_, int vdim_)
474 {
475  if (qspace_ != qspace)
476  {
477  if (own_qspace) { delete qspace; }
478  qspace = qspace_;
479  own_qspace = false;
480  }
481  vdim = (vdim_ < 0) ? vdim : vdim_;
483 }
484 
486  double *qf_data, int vdim_)
487 {
488  if (qspace_ != qspace)
489  {
490  if (own_qspace) { delete qspace; }
491  qspace = qspace_;
492  own_qspace = false;
493  }
494  vdim = (vdim_ < 0) ? vdim : vdim_;
495  NewDataAndSize(qf_data, vdim*qspace->GetSize());
496 }
497 
498 inline void QuadratureFunction::GetElementValues(int idx, Vector &values)
499 {
500  const int s_offset = qspace->element_offsets[idx];
501  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
502  values.NewDataAndSize(data + vdim*s_offset, vdim*sl_size);
503 }
504 
506 {
507  const int s_offset = qspace->element_offsets[idx];
508  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
509  values.Reset(data + vdim*s_offset, vdim, sl_size);
510 }
511 
512 } // namespace mfem
513 
514 #endif
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:2277
int GetVDim() const
Get the vector dimension.
Definition: gridfunc.hpp:398
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.hpp:325
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8491
QuadratureFunction(QuadratureSpace *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace.
Definition: gridfunc.hpp:352
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:94
double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:249
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:242
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:453
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1169
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains)
Definition: gridfunc.cpp:2428
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:368
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:310
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:2582
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:409
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh)
Definition: gridfunc.cpp:856
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: gridfunc.cpp:244
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1361
bool own_qspace
QuadratureSpace ownership flag.
Definition: gridfunc.hpp:341
void GetBdrValuesFrom(GridFunction &)
Definition: gridfunc.cpp:619
void GetGradient(ElementTransformation &tr, Vector &grad)
Definition: gridfunc.cpp:964
int VectorDim() const
Definition: gridfunc.cpp:259
virtual ~QuadratureFunction()
Definition: gridfunc.hpp:368
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
Delta function coefficient.
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
Definition: gridfunc.cpp:1085
int vdim
Vector dimension.
Definition: gridfunc.hpp:340
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
Definition: gridfunc.hpp:473
bool OwnsSpace()
Get the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:405
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:312
QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace, using the external data...
Definition: gridfunc.hpp:360
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:725
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2185
const FiniteElementSpace * FESpace() const
Definition: gridfunc.hpp:294
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: gridfunc.cpp:2572
double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:232
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:796
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1475
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:197
double GetDivergence(ElementTransformation &tr)
Definition: gridfunc.cpp:875
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:217
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:67
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:656
void SetOwnsSpace(bool own)
Set the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:408
double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:1693
GridFunction(const GridFunction &orig)
Copy constructor.
Definition: gridfunc.hpp:66
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:371
QuadratureFunction()
Create an empty QuadratureFunction.
Definition: gridfunc.hpp:346
void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:179
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition: gridfunc.hpp:401
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:280
void GetCurl(ElementTransformation &tr, Vector &curl)
Definition: gridfunc.cpp:911
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:402
FiniteElementCollection * OwnFEC()
Definition: gridfunc.hpp:93
virtual ~ExtrudeCoefficient()
Definition: gridfunc.hpp:463
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:293
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:412
double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2004
const IntegrationRule & GetElementIntRule(int idx)
Get the IntegrationRule associated with mesh element idx.
Definition: gridfunc.hpp:411
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:1052
double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:257
QuadratureSpace * qspace
Associated QuadratureSpace.
Definition: gridfunc.hpp:339
FiniteElementCollection * fec
Used when the grid function is read from a file.
Definition: gridfunc.hpp:34
void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:144
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:2515
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:539
Class for integration point with weight.
Definition: intrules.hpp:25
void SaveSTL(std::ostream &out, int TimesToRefine=1)
Definition: gridfunc.cpp:2297
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:435
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: gridfunc.hpp:498
void GetElementAverages(GridFunction &avgs)
Definition: gridfunc.cpp:1022
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1546
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:485
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1232
void SetSpace(FiniteElementSpace *f)
Definition: gridfunc.cpp:171
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:698
double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:1898
ExtrudeCoefficient(Mesh *m, Coefficient &s, int _n)
Definition: gridfunc.hpp:460
GridFunction & operator=(double value)
Redefine &#39;=&#39; for GridFunction = constant.
Definition: gridfunc.cpp:2160
Vector data type.
Definition: vector.hpp:36
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad)
Definition: gridfunc.cpp:983
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:380
double * data
Definition: vector.hpp:41
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:370
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Definition: gridfunc.cpp:1008
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: gridfunc.cpp:295
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:336
void GetValuesFrom(GridFunction &)
Definition: gridfunc.cpp:582
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
Definition: gridfunc.cpp:2412
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
Definition: gridfunc.cpp:2200
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Definition: gridfunc.hpp:197
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:353
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:756