MFEM  v4.4.0
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-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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 #ifdef MFEM_USE_ADIOS2
20 #include "../general/adios2stream.hpp"
21 #endif
22 #include <limits>
23 #include <ostream>
24 #include <string>
25 
26 namespace mfem
27 {
28 
29 /// Class for grid function - Vector with associated FE space.
30 class GridFunction : public Vector
31 {
32 protected:
33  /// FE space on which the grid function lives. Owned if #fec is not NULL.
35 
36  /** @brief Used when the grid function is read from a file. It can also be
37  set explicitly, see MakeOwner().
38 
39  If not NULL, this pointer is owned by the GridFunction. */
41 
42  long fes_sequence; // see FiniteElementSpace::sequence, Mesh::sequence
43 
44  /** Optional, internal true-dof vector: if the FiniteElementSpace #fes has a
45  non-trivial (i.e. not NULL) prolongation operator, this Vector may hold
46  associated true-dof values - either owned or external. */
48 
49  void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[]);
50 
52 
53  // Project the delta coefficient without scaling and return the (local)
54  // integral of the projection.
55  void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff,
56  double &integral);
57 
58  // Sum fluxes to vertices and count element contributions
60  GridFunction &flux,
61  Array<int>& counts,
62  bool wcoef,
63  int subdomain);
64 
65  /** Project a discontinuous vector coefficient in a continuous space and
66  return in dof_attr the maximal attribute of the elements containing each
67  degree of freedom. */
68  void ProjectDiscCoefficient(VectorCoefficient &coeff, Array<int> &dof_attr);
69 
70  /// Loading helper.
71  void LegacyNCReorder();
72 
73  void Destroy();
74 
75 public:
76 
77  GridFunction() { fes = NULL; fec = NULL; fes_sequence = 0; UseDevice(true); }
78 
79  /// Copy constructor. The internal true-dof vector #t_vec is not copied.
81  : Vector(orig), fes(orig.fes), fec(NULL), fes_sequence(orig.fes_sequence)
82  { UseDevice(true); }
83 
84  /// Construct a GridFunction associated with the FiniteElementSpace @a *f.
86  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
87 
88  /// Construct a GridFunction using previously allocated array @a data.
89  /** The GridFunction does not assume ownership of @a data which is assumed to
90  be of size at least `f->GetVSize()`. Similar to the Vector constructor
91  for externally allocated array, the pointer @a data can be NULL. The data
92  array can be replaced later using the method SetData().
93  */
95  : Vector(data, f->GetVSize())
96  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
97 
98  /** @brief Construct a GridFunction using previously allocated Vector @a base
99  starting at the given offset, @a base_offset. */
100  GridFunction(FiniteElementSpace *f, Vector &base, int base_offset = 0)
101  : Vector(base, base_offset, f->GetVSize())
102  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
103 
104  /// Construct a GridFunction on the given Mesh, using the data from @a input.
105  /** The content of @a input should be in the format created by the method
106  Save(). The reconstructed FiniteElementSpace and FiniteElementCollection
107  are owned by the GridFunction. */
108  GridFunction(Mesh *m, std::istream &input);
109 
110  GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces);
111 
112  /// Copy assignment. Only the data of the base class Vector is copied.
113  /** It is assumed that this object and @a rhs use FiniteElementSpace%s that
114  have the same size.
115 
116  @note Defining this method overwrites the implicitly defined copy
117  assignment operator. */
119  { return operator=((const Vector &)rhs); }
120 
121  /// Make the GridFunction the owner of #fec and #fes.
122  /** If the new FiniteElementCollection, @a fec_, is NULL, ownership of #fec
123  and #fes is taken away. */
124  void MakeOwner(FiniteElementCollection *fec_) { fec = fec_; }
125 
127 
128  int VectorDim() const;
129  int CurlDim() const;
130 
131  /// Read only access to the (optional) internal true-dof Vector.
132  /** Note that the returned Vector may be empty, if not previously allocated
133  or set. */
134  const Vector &GetTrueVector() const { return t_vec; }
135  /// Read and write access to the (optional) internal true-dof Vector.
136  /** Note that the returned Vector may be empty, if not previously allocated
137  or set. */
138  Vector &GetTrueVector() { return t_vec; }
139 
140  /// Extract the true-dofs from the GridFunction.
141  void GetTrueDofs(Vector &tv) const;
142 
143  /// Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
145 
146  /// Set the GridFunction from the given true-dof vector.
147  virtual void SetFromTrueDofs(const Vector &tv);
148 
149  /// Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
151 
152  /// Returns the values in the vertices of i'th element for dimension vdim.
153  void GetNodalValues(int i, Array<double> &nval, int vdim = 1) const;
154 
155  /** @name Element index Get Value Methods
156 
157  These methods take an element index and return the interpolated value of
158  the field at a given reference point within the element.
159 
160  @warning These methods retrieve and use the ElementTransformation object
161  from the mfem::Mesh. This can alter the state of the element
162  transformation object and can also lead to unexpected results when the
163  ElementTransformation object is already in use such as when these methods
164  are called from within an integration loop. Consider using
165  GetValue(ElementTransformation &T, ...) instead.
166  */
167  ///@{
168  /** Return a scalar value from within the given element. */
169  virtual double GetValue(int i, const IntegrationPoint &ip,
170  int vdim = 1) const;
171 
172  /** Return a vector value from within the given element. */
173  virtual void GetVectorValue(int i, const IntegrationPoint &ip,
174  Vector &val) const;
175  ///@}
176 
177  /** @name Element Index Get Values Methods
178 
179  These are convenience methods for repeatedly calling GetValue for
180  multiple points within a given element. The GetValues methods are
181  optimized and should perform better than repeatedly calling GetValue. The
182  GetVectorValues method simply calls GetVectorValue repeatedly.
183 
184  @warning These methods retrieve and use the ElementTransformation object
185  from the mfem::Mesh. This can alter the state of the element
186  transformation object and can also lead to unexpected results when the
187  ElementTransformation object is already in use such as when these methods
188  are called from within an integration loop. Consider using
189  GetValues(ElementTransformation &T, ...) instead.
190  */
191  ///@{
192  /** Compute a collection of scalar values from within the element indicated
193  by the index i. */
194  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
195  int vdim = 1) const;
196 
197  /** Compute a collection of vector values from within the element indicated
198  by the index i. */
199  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
200  DenseMatrix &tr, int vdim = 1) const;
201 
202  void GetVectorValues(int i, const IntegrationRule &ir,
203  DenseMatrix &vals, DenseMatrix &tr) const;
204  ///@}
205 
206  /** @name ElementTransformation Get Value Methods
207 
208  These member functions are designed for use within
209  GridFunctionCoefficient objects. These can be used with
210  ElementTransformation objects coming from either
211  Mesh::GetElementTransformation() or Mesh::GetBdrElementTransformation().
212 
213  @note These methods do not reset the ElementTransformation object so they
214  should be safe to use within integration loops or other contexts where
215  the ElementTransformation is already in use.
216  */
217  ///@{
218  /** Return a scalar value from within the element indicated by the
219  ElementTransformation Object. */
220  virtual double GetValue(ElementTransformation &T, const IntegrationPoint &ip,
221  int comp = 0, Vector *tr = NULL) const;
222 
223  /** Return a vector value from within the element indicated by the
224  ElementTransformation Object. */
225  virtual void GetVectorValue(ElementTransformation &T,
226  const IntegrationPoint &ip,
227  Vector &val, Vector *tr = NULL) const;
228  ///@}
229 
230  /** @name ElementTransformation Get Values Methods
231 
232  These are convenience methods for repeatedly calling GetValue for
233  multiple points within a given element. They work by calling either the
234  ElementTransformation or FaceElementTransformations versions described
235  above. Consequently, these methods should not be expected to run faster
236  than calling the above methods in an external loop.
237 
238  @note These methods do not reset the ElementTransformation object so they
239  should be safe to use within integration loops or other contexts where
240  the ElementTransformation is already in use.
241 
242  @note These methods can also be used with FaceElementTransformations
243  objects.
244  */
245  ///@{
246  /** Compute a collection of scalar values from within the element indicated
247  by the ElementTransformation object. */
249  Vector &vals, int comp = 0, DenseMatrix *tr = NULL) const;
250 
251  /** Compute a collection of vector values from within the element indicated
252  by the ElementTransformation object. */
254  DenseMatrix &vals, DenseMatrix *tr = NULL) const;
255  ///@}
256 
257  /** @name Face Index Get Values Methods
258 
259  These methods are designed to work with Discontinuous Galerkin basis
260  functions. They compute field values on the interface between elements,
261  or on boundary elements, by interpolating the field in a neighboring
262  element. The \a side argument indices which neighboring element should be
263  used: 0, 1, or 2 (automatically chosen).
264 
265  @warning These methods retrieve and use the FaceElementTransformations
266  object from the mfem::Mesh. This can alter the state of the face element
267  transformations object and can also lead to unexpected results when the
268  FaceElementTransformations object is already in use such as when these
269  methods are called from within an integration loop. Consider using
270  GetValues(ElementTransformation &T, ...) instead.
271  */
272  ///@{
273  /** Compute a collection of scalar values from within the face
274  indicated by the index i. */
275  int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals,
276  DenseMatrix &tr, int vdim = 1) const;
277 
278  /** Compute a collection of vector values from within the face
279  indicated by the index i. */
280  int GetFaceVectorValues(int i, int side, const IntegrationRule &ir,
281  DenseMatrix &vals, DenseMatrix &tr) const;
282  ///@}
283 
284  void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps,
285  int vdim = 1) const;
286 
287  void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps,
288  DenseMatrix &tr, int vdim = 1) const;
289 
290  void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess,
291  int vdim = 1) const;
292 
293  void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess,
294  DenseMatrix &tr, int vdim = 1) const;
295 
296  void GetValuesFrom(const GridFunction &orig_func);
297 
298  void GetBdrValuesFrom(const GridFunction &orig_func);
299 
300  void GetVectorFieldValues(int i, const IntegrationRule &ir,
301  DenseMatrix &vals,
302  DenseMatrix &tr, int comp = 0) const;
303 
304  /// For a vector grid function, makes sure that the ordering is byNODES.
305  void ReorderByNodes();
306 
307  /// Return the values as a vector on mesh vertices for dimension vdim.
308  void GetNodalValues(Vector &nval, int vdim = 1) const;
309 
310  void GetVectorFieldNodalValues(Vector &val, int comp) const;
311 
312  void ProjectVectorFieldOn(GridFunction &vec_field, int comp = 0);
313 
314  /** @brief Compute a certain derivative of a function's component.
315  Derivatives of the function are computed at the DOF locations of @a der,
316  and averaged over overlapping DOFs. Thus this function projects the
317  derivative to the FiniteElementSpace of @a der.
318  @param[in] comp Index of the function's component to be differentiated.
319  The index is 1-based, i.e., use 1 for scalar functions.
320  @param[in] der_comp Use 0/1/2 for derivatives in x/y/z directions.
321  @param[out] der The resulting derivative (scalar function). The
322  FiniteElementSpace of this function must be set
323  before the call. */
324  void GetDerivative(int comp, int der_comp, GridFunction &der);
325 
326  double GetDivergence(ElementTransformation &tr) const;
327 
328  void GetCurl(ElementTransformation &tr, Vector &curl) const;
329 
330  void GetGradient(ElementTransformation &tr, Vector &grad) const;
331 
333  DenseMatrix &grad) const;
334 
335  void GetGradients(const int elem, const IntegrationRule &ir,
336  DenseMatrix &grad) const
337  { GetGradients(*fes->GetElementTransformation(elem), ir, grad); }
338 
339  void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const;
340 
341  /** Compute \f$ (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) \f$,
342  where \f$ \psi_i \f$ are the basis functions for the FE space of avgs.
343  Both FE spaces should be scalar and on the same mesh. */
344  void GetElementAverages(GridFunction &avgs) const;
345 
346  /** Sets the output vector @a dof_vals to the values of the degrees of
347  freedom of element @a el. */
348  virtual void GetElementDofValues(int el, Vector &dof_vals) const;
349 
350  /** Impose the given bounds on the function's DOFs while preserving its local
351  * integral (described in terms of the given weights) on the i'th element
352  * through SLBPQ optimization.
353  * Intended to be used for discontinuous FE functions. */
354  void ImposeBounds(int i, const Vector &weights,
355  const Vector &lo_, const Vector &hi_);
356  void ImposeBounds(int i, const Vector &weights,
357  double min_ = 0.0, double max_ = infinity());
358 
359  /** On a non-conforming mesh, make sure the function lies in the conforming
360  space by multiplying with R and then with P, the conforming restriction
361  and prolongation matrices of the space, respectively. */
362  void RestrictConforming();
363 
364  /** @brief Project the @a src GridFunction to @a this GridFunction, both of
365  which must be on the same mesh. */
366  /** The current implementation assumes that all elements use the same
367  projection matrix. */
368  void ProjectGridFunction(const GridFunction &src);
369 
370  /** @brief Project @a coeff Coefficient to @a this GridFunction. The
371  projection computation depends on the choice of the FiniteElementSpace
372  #fes. Note that this is usually interpolation at the degrees of freedom
373  in each element (not L2 projection). */
374  virtual void ProjectCoefficient(Coefficient &coeff);
375 
376  /** @brief Project @a coeff Coefficient to @a this GridFunction, using one
377  element for each degree of freedom in @a dofs and nodal interpolation on
378  that element. */
379  void ProjectCoefficient(Coefficient &coeff, Array<int> &dofs, int vd = 0);
380 
381  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction. The
382  projection computation depends on the choice of the FiniteElementSpace
383  #fes. Note that this is usually interpolation at the degrees of freedom
384  in each element (not L2 projection).*/
386 
387  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction, using
388  one element for each degree of freedom in @a dofs and nodal interpolation
389  on that element. */
390  void ProjectCoefficient(VectorCoefficient &vcoeff, Array<int> &dofs);
391 
392  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction, only
393  projecting onto elements with the given @a attribute */
394  void ProjectCoefficient(VectorCoefficient &vcoeff, int attribute);
395 
396  /** @brief Analogous to the version with argument @a vcoeff VectorCoefficient
397  but using an array of scalar coefficients for each component. */
398  void ProjectCoefficient(Coefficient *coeff[]);
399 
400  /** @brief Project a discontinuous vector coefficient as a grid function on
401  a continuous finite element space. The values in shared dofs are
402  determined from the element with maximal attribute. */
403  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
404 
406  /** @brief Projects a discontinuous coefficient so that the values in shared
407  vdofs are computed by taking an average of the possible values. */
408  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
409  /** @brief Projects a discontinuous _vector_ coefficient so that the values
410  in shared vdofs are computed by taking an average of the possible values.
411  */
412  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff, AvgType type);
413 
414 protected:
415  /** @brief Accumulates (depending on @a type) the values of @a coeff at all
416  shared vdofs and counts in how many zones each vdof appears. */
417  void AccumulateAndCountZones(Coefficient &coeff, AvgType type,
418  Array<int> &zones_per_vdof);
419 
420  /** @brief Accumulates (depending on @a type) the values of @a vcoeff at all
421  shared vdofs and counts in how many zones each vdof appears. */
423  Array<int> &zones_per_vdof);
424 
425  /** @brief Used for the serial and parallel implementations of the
426  GetDerivative() method; see its documentation. */
427  void AccumulateAndCountDerivativeValues(int comp, int der_comp,
428  GridFunction &der,
429  Array<int> &zones_per_dof);
430 
432  VectorCoefficient *vcoeff, Array<int> &attr,
433  Array<int> &values_counter);
434 
436  Array<int> &bdr_attr,
437  Array<int> &values_counter);
438 
439  // Complete the computation of averages; called e.g. after
440  // AccumulateAndCountZones().
441  void ComputeMeans(AvgType type, Array<int> &zones_per_vdof);
442 
443 public:
444  /** @brief Project a Coefficient on the GridFunction, modifying only DOFs on
445  the boundary associated with the boundary attributes marked in the
446  @a attr array. */
448  {
449  Coefficient *coeff_p = &coeff;
450  ProjectBdrCoefficient(&coeff_p, attr);
451  }
452 
453  /** @brief Project a VectorCoefficient on the GridFunction, modifying only
454  DOFs on the boundary associated with the boundary attributes marked in
455  the @a attr array. */
456  virtual void ProjectBdrCoefficient(VectorCoefficient &vcoeff,
457  Array<int> &attr);
458 
459  /** @brief Project a set of Coefficient%s on the components of the
460  GridFunction, modifying only DOFs on the boundary associated with the
461  boundary attributed marked in the @a attr array. */
462  /** If a Coefficient pointer in the array @a coeff is NULL, that component
463  will not be touched. */
464  virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
465 
466  /** Project the normal component of the given VectorCoefficient on
467  the boundary. Only boundary attributes that are marked in
468  'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
470  Array<int> &bdr_attr);
471 
472  /** @brief Project the tangential components of the given VectorCoefficient
473  on the boundary. Only boundary attributes that are marked in @a bdr_attr
474  are projected. Assumes ND-type VectorFE GridFunction. */
475  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
476  Array<int> &bdr_attr);
477 
478 
479  virtual double ComputeL2Error(Coefficient &exsol,
480  const IntegrationRule *irs[] = NULL) const
481  { return ComputeLpError(2.0, exsol, NULL, irs); }
482 
483  virtual double ComputeL2Error(Coefficient *exsol[],
484  const IntegrationRule *irs[] = NULL) const;
485 
486  virtual double ComputeL2Error(VectorCoefficient &exsol,
487  const IntegrationRule *irs[] = NULL,
488  Array<int> *elems = NULL) const;
489 
490  /// Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements
491  virtual double ComputeGradError(VectorCoefficient *exgrad,
492  const IntegrationRule *irs[] = NULL) const;
493 
494  /// Returns ||curl u_ex - curl u_h||_L2 for ND elements
495  virtual double ComputeCurlError(VectorCoefficient *excurl,
496  const IntegrationRule *irs[] = NULL) const;
497 
498  /// Returns ||div u_ex - div u_h||_L2 for RT elements
499  virtual double ComputeDivError(Coefficient *exdiv,
500  const IntegrationRule *irs[] = NULL) const;
501 
502  /// Returns the Face Jumps error for L2 elements. The error can be weighted
503  /// by a constant nu, by nu/h, or nu*p^2/h, depending on the value of
504  /// @a jump_scaling.
505  virtual double ComputeDGFaceJumpError(Coefficient *exsol,
506  Coefficient *ell_coeff,
507  class JumpScaling jump_scaling,
508  const IntegrationRule *irs[] = NULL)
509  const;
510 
511  /// Returns the Face Jumps error for L2 elements, with 1/h scaling.
512  MFEM_DEPRECATED
513  double ComputeDGFaceJumpError(Coefficient *exsol,
514  Coefficient *ell_coeff,
515  double Nu,
516  const IntegrationRule *irs[] = NULL) const;
517 
518  /** This method is kept for backward compatibility.
519 
520  Returns either the H1-seminorm, or the DG face jumps error, or both
521  depending on norm_type = 1, 2, 3. Additional arguments for the DG face
522  jumps norm: ell_coeff: mesh-depended coefficient (weight) Nu: scalar
523  constant weight */
524  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
525  Coefficient *ell_coef, double Nu,
526  int norm_type) const;
527 
528  /// Returns the error measured in H1-norm for H1 elements or in "broken"
529  /// H1-norm for L2 elements
530  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
531  const IntegrationRule *irs[] = NULL) const;
532 
533  /// Returns the error measured in H(div)-norm for RT elements
534  virtual double ComputeHDivError(VectorCoefficient *exsol,
535  Coefficient *exdiv,
536  const IntegrationRule *irs[] = NULL) const;
537 
538  /// Returns the error measured in H(curl)-norm for ND elements
539  virtual double ComputeHCurlError(VectorCoefficient *exsol,
540  VectorCoefficient *excurl,
541  const IntegrationRule *irs[] = NULL) const;
542 
543  virtual double ComputeMaxError(Coefficient &exsol,
544  const IntegrationRule *irs[] = NULL) const
545  {
546  return ComputeLpError(infinity(), exsol, NULL, irs);
547  }
548 
549  virtual double ComputeMaxError(Coefficient *exsol[],
550  const IntegrationRule *irs[] = NULL) const;
551 
552  virtual double ComputeMaxError(VectorCoefficient &exsol,
553  const IntegrationRule *irs[] = NULL) const
554  {
555  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
556  }
557 
558  virtual double ComputeL1Error(Coefficient &exsol,
559  const IntegrationRule *irs[] = NULL) const
560  { return ComputeLpError(1.0, exsol, NULL, irs); }
561 
562  virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
563  int norm_type, Array<int> *elems = NULL,
564  const IntegrationRule *irs[] = NULL) const;
565 
566  virtual double ComputeL1Error(VectorCoefficient &exsol,
567  const IntegrationRule *irs[] = NULL) const
568  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
569 
570  virtual double ComputeLpError(const double p, Coefficient &exsol,
571  Coefficient *weight = NULL,
572  const IntegrationRule *irs[] = NULL) const;
573 
574  /** Compute the Lp error in each element of the mesh and store the results in
575  the Vector @a error. The result should be of length number of elements,
576  for example an L2 GridFunction of order zero using map type VALUE. */
577  virtual void ComputeElementLpErrors(const double p, Coefficient &exsol,
578  Vector &error,
579  Coefficient *weight = NULL,
580  const IntegrationRule *irs[] = NULL
581  ) const;
582 
583  virtual void ComputeElementL1Errors(Coefficient &exsol,
584  Vector &error,
585  const IntegrationRule *irs[] = NULL
586  ) const
587  { ComputeElementLpErrors(1.0, exsol, error, NULL, irs); }
588 
589  virtual void ComputeElementL2Errors(Coefficient &exsol,
590  Vector &error,
591  const IntegrationRule *irs[] = NULL
592  ) const
593  { ComputeElementLpErrors(2.0, exsol, error, NULL, irs); }
594 
596  Vector &error,
597  const IntegrationRule *irs[] = NULL
598  ) const
599  { ComputeElementLpErrors(infinity(), exsol, error, NULL, irs); }
600 
601  /** When given a vector weight, compute the pointwise (scalar) error as the
602  dot product of the vector error with the vector weight. Otherwise, the
603  scalar error is the l_2 norm of the vector error. */
604  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
605  Coefficient *weight = NULL,
606  VectorCoefficient *v_weight = NULL,
607  const IntegrationRule *irs[] = NULL) const;
608 
609  /** Compute the Lp error in each element of the mesh and store the results in
610  the Vector @ error. The result should be of length number of elements,
611  for example an L2 GridFunction of order zero using map type VALUE. */
612  virtual void ComputeElementLpErrors(const double p, VectorCoefficient &exsol,
613  Vector &error,
614  Coefficient *weight = NULL,
615  VectorCoefficient *v_weight = NULL,
616  const IntegrationRule *irs[] = NULL
617  ) const;
618 
620  Vector &error,
621  const IntegrationRule *irs[] = NULL
622  ) const
623  { ComputeElementLpErrors(1.0, exsol, error, NULL, NULL, irs); }
624 
626  Vector &error,
627  const IntegrationRule *irs[] = NULL
628  ) const
629  { ComputeElementLpErrors(2.0, exsol, error, NULL, NULL, irs); }
630 
632  Vector &error,
633  const IntegrationRule *irs[] = NULL
634  ) const
635  { ComputeElementLpErrors(infinity(), exsol, error, NULL, NULL, irs); }
636 
637  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
638  GridFunction &flux,
639  bool wcoef = true, int subdomain = -1);
640 
641  /// Redefine '=' for GridFunction = constant.
642  GridFunction &operator=(double value);
643 
644  /// Copy the data from @a v.
645  /** The size of @a v must be equal to the size of the associated
646  FiniteElementSpace #fes. */
647  GridFunction &operator=(const Vector &v);
648 
649  /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
650  virtual void Update();
651 
653  const FiniteElementSpace *FESpace() const { return fes; }
654 
655  /// Associate a new FiniteElementSpace with the GridFunction.
656  /** The GridFunction is resized using the SetSize() method. */
657  virtual void SetSpace(FiniteElementSpace *f);
658 
659  using Vector::MakeRef;
660 
661  /** @brief Make the GridFunction reference external data on a new
662  FiniteElementSpace. */
663  /** This method changes the FiniteElementSpace associated with the
664  GridFunction and sets the pointer @a v as external data in the
665  GridFunction. */
666  virtual void MakeRef(FiniteElementSpace *f, double *v);
667 
668  /** @brief Make the GridFunction reference external data on a new
669  FiniteElementSpace. */
670  /** This method changes the FiniteElementSpace associated with the
671  GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
672  as external data in the GridFunction.
673  @note This version of the method will also perform bounds checks when
674  the build option MFEM_DEBUG is enabled. */
675  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
676 
677  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
678  GridFunction. */
679  /** - If the prolongation matrix of @a f is trivial (i.e. its method
680  FiniteElementSpace::GetProlongationMatrix() returns NULL), then the
681  method MakeRef() is called with the same arguments.
682  - Otherwise, the method SetSpace() is called with argument @a f.
683  - The internal true-dof vector is set to reference @a tv. */
684  void MakeTRef(FiniteElementSpace *f, double *tv);
685 
686  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
687  GridFunction. */
688  /** - If the prolongation matrix of @a f is trivial (i.e. its method
689  FiniteElementSpace::GetProlongationMatrix() returns NULL), this method
690  calls MakeRef() with the same arguments.
691  - Otherwise, this method calls SetSpace() with argument @a f.
692  - The internal true-dof vector is set to reference the sub-vector of
693  @a tv starting at the offset @a tv_offset. */
694  void MakeTRef(FiniteElementSpace *f, Vector &tv, int tv_offset);
695 
696  /// Save the GridFunction to an output stream.
697  virtual void Save(std::ostream &out) const;
698 
699  /// Save the GridFunction to a file. The given @a precision will be used for
700  /// ASCII output.
701  virtual void Save(const char *fname, int precision=16) const;
702 
703 #ifdef MFEM_USE_ADIOS2
704  /// Save the GridFunction to a binary output stream using adios2 bp format.
705  virtual void Save(adios2stream &out, const std::string& variable_name,
708 #endif
709 
710  /** @brief Write the GridFunction in VTK format. Note that Mesh::PrintVTK
711  must be called first. The parameter ref > 0 must match the one used in
712  Mesh::PrintVTK. */
713  void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
714 
715  /** @brief Write the GridFunction in STL format. Note that the mesh dimension
716  must be 2 and that quad elements will be broken into two triangles.*/
717  void SaveSTL(std::ostream &out, int TimesToRefine = 1);
718 
719  /// Destroys grid function.
720  virtual ~GridFunction() { Destroy(); }
721 };
722 
723 
724 /** Overload operator<< for std::ostream and GridFunction; valid also for the
725  derived class ParGridFunction */
726 std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
727 
728 /// Class used to specify how the jump terms in
729 /// GridFunction::ComputeDGFaceJumpError are scaled.
731 {
732 public:
734  {
738  };
739 private:
740  double nu;
741  JumpScalingType type;
742 public:
743  JumpScaling(double nu_=1.0, JumpScalingType type_=CONSTANT)
744  : nu(nu_), type(type_) { }
745  double Eval(double h, int p) const
746  {
747  double val = nu;
748  if (type != CONSTANT) { val /= h; }
749  if (type == P_SQUARED_OVER_H) { val *= p*p; }
750  return val;
751  }
752 };
753 
754 
755 /** @brief Class representing a function through its values (scalar or vector)
756  at quadrature points. */
758 {
759 protected:
760  QuadratureSpace *qspace; ///< Associated QuadratureSpace
761  int vdim; ///< Vector dimension
762  bool own_qspace; ///< QuadratureSpace ownership flag
763 
764 public:
765  /// Create an empty QuadratureFunction.
766  /** The object can be initialized later using the SetSpace() methods. */
768  : qspace(NULL), vdim(0), own_qspace(false) { }
769 
770  /** @brief Copy constructor. The QuadratureSpace ownership flag, #own_qspace,
771  in the new object is set to false. */
773  : Vector(orig),
774  qspace(orig.qspace), vdim(orig.vdim), own_qspace(false) { }
775 
776  /// Create a QuadratureFunction based on the given QuadratureSpace.
777  /** The QuadratureFunction does not assume ownership of the QuadratureSpace.
778  @note The Vector data is not initialized. */
779  QuadratureFunction(QuadratureSpace *qspace_, int vdim_ = 1)
780  : Vector(vdim_*qspace_->GetSize()),
781  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
782 
783  /** @brief Create a QuadratureFunction based on the given QuadratureSpace,
784  using the external data, @a qf_data. */
785  /** The QuadratureFunction does not assume ownership of neither the
786  QuadratureSpace nor the external data. */
787  QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_ = 1)
788  : Vector(qf_data, vdim_*qspace_->GetSize()),
789  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
790 
791  /// Read a QuadratureFunction from the stream @a in.
792  /** The QuadratureFunction assumes ownership of the read QuadratureSpace. */
793  QuadratureFunction(Mesh *mesh, std::istream &in);
794 
795  virtual ~QuadratureFunction() { if (own_qspace) { delete qspace; } }
796 
797  /// Get the associated QuadratureSpace.
798  QuadratureSpace *GetSpace() const { return qspace; }
799 
800  /// Change the QuadratureSpace and optionally the vector dimension.
801  /** If the new QuadratureSpace is different from the current one, the
802  QuadratureFunction will not assume ownership of the new space; otherwise,
803  the ownership flag remains the same.
804 
805  If the new vector dimension @a vdim_ < 0, the vector dimension remains
806  the same.
807 
808  The data size is updated by calling Vector::SetSize(). */
809  inline void SetSpace(QuadratureSpace *qspace_, int vdim_ = -1);
810 
811  /** @brief Change the QuadratureSpace, the data array, and optionally the
812  vector dimension. */
813  /** If the new QuadratureSpace is different from the current one, the
814  QuadratureFunction will not assume ownership of the new space; otherwise,
815  the ownership flag remains the same.
816 
817  If the new vector dimension @a vdim_ < 0, the vector dimension remains
818  the same.
819 
820  The data array is replaced by calling Vector::NewDataAndSize(). */
821  inline void SetSpace(QuadratureSpace *qspace_, double *qf_data,
822  int vdim_ = -1);
823 
824  /// Get the vector dimension.
825  int GetVDim() const { return vdim; }
826 
827  /// Set the vector dimension, updating the size by calling Vector::SetSize().
828  void SetVDim(int vdim_)
829  { vdim = vdim_; SetSize(vdim*qspace->GetSize()); }
830 
831  /// Get the QuadratureSpace ownership flag.
832  bool OwnsSpace() { return own_qspace; }
833 
834  /// Set the QuadratureSpace ownership flag.
835  void SetOwnsSpace(bool own) { own_qspace = own; }
836 
837  /// Redefine '=' for QuadratureFunction = constant.
838  QuadratureFunction &operator=(double value);
839 
840  /// Copy the data from @a v.
841  /** The size of @a v must be equal to the size of the associated
842  QuadratureSpace #qspace times the QuadratureFunction dimension
843  i.e. QuadratureFunction::Size(). */
845 
846  /// Copy assignment. Only the data of the base class Vector is copied.
847  /** The QuadratureFunctions @a v and @a *this must have QuadratureSpaces with
848  the same size.
849 
850  @note Defining this method overwrites the implicitly defined copy
851  assignment operator. */
853 
854  /// Get the IntegrationRule associated with mesh element @a idx.
855  const IntegrationRule &GetElementIntRule(int idx) const
856  { return qspace->GetElementIntRule(idx); }
857 
858  /// Return all values associated with mesh element @a idx in a Vector.
859  /** The result is stored in the Vector @a values as a reference to the
860  global values.
861 
862  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
863  `i`-th vector component at the `j`-th quadrature point.
864  */
865  inline void GetElementValues(int idx, Vector &values);
866 
867  /// Return all values associated with mesh element @a idx in a Vector.
868  /** The result is stored in the Vector @a values as a copy of the
869  global values.
870 
871  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
872  `i`-th vector component at the `j`-th quadrature point.
873  */
874  inline void GetElementValues(int idx, Vector &values) const;
875 
876  /// Return the quadrature function values at an integration point.
877  /** The result is stored in the Vector @a values as a reference to the
878  global values. */
879  inline void GetElementValues(int idx, const int ip_num, Vector &values);
880 
881  /// Return the quadrature function values at an integration point.
882  /** The result is stored in the Vector @a values as a copy to the
883  global values. */
884  inline void GetElementValues(int idx, const int ip_num, Vector &values) const;
885 
886  /// Return all values associated with mesh element @a idx in a DenseMatrix.
887  /** The result is stored in the DenseMatrix @a values as a reference to the
888  global values.
889 
890  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
891  `i`-th vector component at the `j`-th quadrature point.
892  */
893  inline void GetElementValues(int idx, DenseMatrix &values);
894 
895  /// Return all values associated with mesh element @a idx in a const DenseMatrix.
896  /** The result is stored in the DenseMatrix @a values as a copy of the
897  global values.
898 
899  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
900  `i`-th vector component at the `j`-th quadrature point.
901  */
902  inline void GetElementValues(int idx, DenseMatrix &values) const;
903 
904  /// Write the QuadratureFunction to the stream @a out.
905  void Save(std::ostream &out) const;
906 
907  /// @brief Write the QuadratureFunction to @a out in VTU (ParaView) format.
908  ///
909  /// The data will be uncompressed if @a compression_level is zero, or if the
910  /// format is VTKFormat::ASCII. Otherwise, zlib compression will be used for
911  /// binary data.
912  void SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII,
913  int compression_level=0) const;
914 
915  /// @brief Save the QuadratureFunction to a VTU (ParaView) file.
916  ///
917  /// The extension ".vtu" will be appended to @a filename.
918  /// @sa SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII,
919  /// int compression_level=0)
920  void SaveVTU(const std::string &filename, VTKFormat format=VTKFormat::ASCII,
921  int compression_level=0) const;
922 };
923 
924 /// Overload operator<< for std::ostream and QuadratureFunction.
925 std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
926 
927 
928 double ZZErrorEstimator(BilinearFormIntegrator &blfi,
929  GridFunction &u,
930  GridFunction &flux,
931  Vector &error_estimates,
932  Array<int> *aniso_flags = NULL,
933  int with_subdomains = 1,
934  bool with_coeff = false);
935 
936 /// Compute the Lp distance between two grid functions on the given element.
937 double ComputeElementLpDistance(double p, int i,
938  GridFunction& gf1, GridFunction& gf2);
939 
940 
941 /// Class used for extruding scalar GridFunctions
943 {
944 private:
945  int n;
946  Mesh *mesh_in;
947  Coefficient &sol_in;
948 public:
950  : n(n_), mesh_in(m), sol_in(s) { }
951  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
952  virtual ~ExtrudeCoefficient() { }
953 };
954 
955 /// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
956 GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
957  GridFunction *sol, const int ny);
958 
959 
960 // Inline methods
961 
962 inline void QuadratureFunction::SetSpace(QuadratureSpace *qspace_, int vdim_)
963 {
964  if (qspace_ != qspace)
965  {
966  if (own_qspace) { delete qspace; }
967  qspace = qspace_;
968  own_qspace = false;
969  }
970  vdim = (vdim_ < 0) ? vdim : vdim_;
972 }
973 
975  double *qf_data, int vdim_)
976 {
977  if (qspace_ != qspace)
978  {
979  if (own_qspace) { delete qspace; }
980  qspace = qspace_;
981  own_qspace = false;
982  }
983  vdim = (vdim_ < 0) ? vdim : vdim_;
984  NewDataAndSize(qf_data, vdim*qspace->GetSize());
985 }
986 
987 inline void QuadratureFunction::GetElementValues(int idx, Vector &values)
988 {
989  const int s_offset = qspace->element_offsets[idx];
990  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
991  values.NewDataAndSize(data + vdim*s_offset, vdim*sl_size);
992 }
993 
994 inline void QuadratureFunction::GetElementValues(int idx, Vector &values) const
995 {
996  const int s_offset = qspace->element_offsets[idx];
997  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
998  values.SetSize(vdim*sl_size);
999  const double *q = data + vdim*s_offset;
1000  for (int i = 0; i<values.Size(); i++)
1001  {
1002  values(i) = *(q++);
1003  }
1004 }
1005 
1006 inline void QuadratureFunction::GetElementValues(int idx, const int ip_num,
1007  Vector &values)
1008 {
1009  const int s_offset = qspace->element_offsets[idx] * vdim + ip_num * vdim;
1010  values.NewDataAndSize(data + s_offset, vdim);
1011 }
1012 
1013 inline void QuadratureFunction::GetElementValues(int idx, const int ip_num,
1014  Vector &values) const
1015 {
1016  const int s_offset = qspace->element_offsets[idx] * vdim + ip_num * vdim;
1017  values.SetSize(vdim);
1018  const double *q = data + s_offset;
1019  for (int i = 0; i < values.Size(); i++)
1020  {
1021  values(i) = *(q++);
1022  }
1023 }
1024 
1026 {
1027  const int s_offset = qspace->element_offsets[idx];
1028  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
1029  values.Reset(data + vdim*s_offset, vdim, sl_size);
1030 }
1031 
1033  DenseMatrix &values) const
1034 {
1035  const int s_offset = qspace->element_offsets[idx];
1036  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
1037  values.SetSize(vdim, sl_size);
1038  const double *q = data + vdim*s_offset;
1039  for (int j = 0; j<sl_size; j++)
1040  {
1041  for (int i = 0; i<vdim; i++)
1042  {
1043  values(i,j) = *(q++);
1044  }
1045  }
1046 }
1047 
1048 } // namespace mfem
1049 
1050 #endif
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:3739
int GetVDim() const
Get the vector dimension.
Definition: gridfunc.hpp:825
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.hpp:720
QuadratureFunction(QuadratureSpace *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace.
Definition: gridfunc.hpp:779
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void GetElementAverages(GridFunction &avgs) const
Definition: gridfunc.cpp:1800
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:703
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:162
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: gridfunc.hpp:855
Memory< double > data
Definition: vector.hpp:64
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:942
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
Definition: gridfunc.cpp:249
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:150
Base class for vector Coefficients that optionally depend on time and space.
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:479
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:2332
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:976
int CurlDim() const
Definition: gridfunc.cpp:344
virtual void ComputeElementLpErrors(const double p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3356
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: gridfunc.cpp:2885
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:465
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
void RestrictConforming()
Definition: gridfunc.cpp:1978
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:4251
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2569
bool own_qspace
QuadratureSpace ownership flag.
Definition: gridfunc.hpp:762
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:964
int VectorDim() const
Definition: gridfunc.cpp:321
virtual ~QuadratureFunction()
Definition: gridfunc.hpp:795
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, Array< int > &bdr_attr, Array< int > &values_counter)
Definition: gridfunc.cpp:2255
virtual double ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
Definition: gridfunc.cpp:3111
int vdim
Vector dimension.
Definition: gridfunc.hpp:761
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Definition: gridfunc.cpp:1433
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:138
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:124
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:552
ExtrudeCoefficient(Mesh *m, Coefficient &s, int n_)
Definition: gridfunc.hpp:949
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
Definition: gridfunc.hpp:962
Data arrays will be written in ASCII format.
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:2307
bool OwnsSpace()
Get the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:832
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:393
virtual void ComputeElementL2Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:589
QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace, using the external data...
Definition: gridfunc.hpp:787
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:1297
void ImposeBounds(int i, const Vector &weights, const Vector &lo_, const Vector &hi_)
Definition: gridfunc.cpp:1903
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3619
const FiniteElementSpace * FESpace() const
Definition: gridfunc.hpp:653
virtual void GetElementDofValues(int el, Vector &dof_vals) const
Definition: gridfunc.cpp:1840
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Definition: gridfunc.cpp:4241
GridFunction(FiniteElementSpace *f, double *data)
Construct a GridFunction using previously allocated array data.
Definition: gridfunc.hpp:94
void GetDerivative(int comp, int der_comp, GridFunction &der)
Compute a certain derivative of a function&#39;s component. Derivatives of the function are computed at t...
Definition: gridfunc.cpp:1421
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:2663
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:144
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess, int vdim=1) const
Definition: gridfunc.cpp:597
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:558
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
Definition: gridfunc.hpp:335
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:90
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:1232
void SetOwnsSpace(bool own)
Set the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:835
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:3085
GridFunction(const GridFunction &orig)
Copy constructor. The internal true-dof vector t_vec is not copied.
Definition: gridfunc.hpp:80
JumpScaling(double nu_=1.0, JumpScalingType type_=CONSTANT)
Definition: gridfunc.hpp:743
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:219
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:798
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:566
QuadratureFunction()
Create an empty QuadratureFunction.
Definition: gridfunc.hpp:767
VTKFormat
Data array format for VTK and VTU files.
Definition: vtk.hpp:93
void AccumulateAndCountDerivativeValues(int comp, int der_comp, GridFunction &der, Array< int > &zones_per_dof)
Used for the serial and parallel implementations of the GetDerivative() method; see its documentation...
Definition: gridfunc.cpp:1368
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:2017
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:201
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition: gridfunc.hpp:828
void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps, int vdim=1) const
Definition: gridfunc.cpp:558
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition: gridfunc.cpp:364
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:510
FiniteElementCollection * OwnFEC()
Definition: gridfunc.hpp:126
virtual ~ExtrudeCoefficient()
Definition: gridfunc.hpp:952
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: gridfunc.hpp:118
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3291
virtual void ComputeElementL1Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:583
QuadratureFunction & operator=(double value)
Redefine &#39;=&#39; for QuadratureFunction = constant.
Definition: gridfunc.cpp:3936
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:134
void SaveVTU(std::ostream &out, VTKFormat format=VTKFormat::ASCII, int compression_level=0) const
Write the QuadratureFunction to out in VTU (ParaView) format.
Definition: gridfunc.cpp:3969
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:623
virtual double ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
Definition: gridfunc.cpp:2845
virtual void ComputeElementL2Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:625
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
virtual void ComputeElementL1Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:619
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
FiniteElementSpace * fes
FE space on which the grid function lives. Owned if fec is not NULL.
Definition: gridfunc.hpp:34
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
Definition: gridfunc.cpp:1851
QuadratureSpace * qspace
Associated QuadratureSpace.
Definition: gridfunc.hpp:760
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:543
FiniteElementCollection * fec
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner()...
Definition: gridfunc.hpp:40
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1738
double GetDivergence(ElementTransformation &tr) const
Definition: gridfunc.cpp:1456
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:164
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:4184
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: gridfunc.cpp:2926
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:1103
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
Definition: gridfunc.cpp:3120
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains, bool with_coeff)
Definition: gridfunc.cpp:4093
virtual void ComputeElementMaxErrors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:595
Class for integration point with weight.
Definition: intrules.hpp:25
void SaveSTL(std::ostream &out, int TimesToRefine=1)
Write the GridFunction in STL format. Note that the mesh dimension must be 2 and that quad elements w...
Definition: gridfunc.cpp:3759
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:653
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: gridfunc.hpp:987
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Definition: gridfunc.cpp:1709
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
Definition: gridfunc.cpp:2734
void LegacyNCReorder()
Loading helper.
Definition: gridfunc.cpp:3859
QuadratureFunction(const QuadratureFunction &orig)
Copy constructor. The QuadratureSpace ownership flag, own_qspace, in the new object is set to false...
Definition: gridfunc.hpp:772
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:46
GridFunction(FiniteElementSpace *f, Vector &base, int base_offset=0)
Construct a GridFunction using previously allocated Vector base starting at the given offset...
Definition: gridfunc.hpp:100
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2395
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:193
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:1270
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:680
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3185
void GetCurl(ElementTransformation &tr, Vector &curl) const
Definition: gridfunc.cpp:1546
Vector data type.
Definition: vector.hpp:60
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:585
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:120
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: gridfunc.cpp:306
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Definition: gridfunc.hpp:85
RefCoord s[3]
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2962
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:935
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
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:66
long GetSequence() const
Definition: fespace.hpp:898
double Eval(double h, int p) const
Definition: gridfunc.hpp:745
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: gridfunc.cpp:379
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:757
void GetValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:1146
void GetBdrValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:1193
virtual void ComputeElementMaxErrors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:631
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr, Array< int > &values_counter)
Definition: gridfunc.cpp:2113
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
Definition: gridfunc.cpp:3954
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be called first. The parameter ref &gt; 0 must match the one used in Mesh::PrintVTK.
Definition: gridfunc.cpp:3660
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:447
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:438
void GetGradient(ElementTransformation &tr, Vector &grad) const
Definition: gridfunc.cpp:1641
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:1328