MFEM  v4.6.0
Finite element discretization library
gridfunc.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 
51  // Project the delta coefficient without scaling and return the (local)
52  // integral of the projection.
53  void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff,
54  double &integral);
55 
56  // Sum fluxes to vertices and count element contributions
58  GridFunction &flux,
59  Array<int>& counts,
60  bool wcoef,
61  int subdomain);
62 
63  /** Project a discontinuous vector coefficient in a continuous space and
64  return in dof_attr the maximal attribute of the elements containing each
65  degree of freedom. */
66  void ProjectDiscCoefficient(VectorCoefficient &coeff, Array<int> &dof_attr);
67 
68  /// Loading helper.
69  void LegacyNCReorder();
70 
71  void Destroy();
72 
73 public:
74 
75  GridFunction() { fes = NULL; fec = NULL; fes_sequence = 0; UseDevice(true); }
76 
77  /// Copy constructor. The internal true-dof vector #t_vec is not copied.
79  : Vector(orig), fes(orig.fes), fec(NULL), fes_sequence(orig.fes_sequence)
80  { UseDevice(true); }
81 
82  /// Construct a GridFunction associated with the FiniteElementSpace @a *f.
84  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
85 
86  /// Construct a GridFunction using previously allocated array @a data.
87  /** The GridFunction does not assume ownership of @a data which is assumed to
88  be of size at least `f->GetVSize()`. Similar to the Vector constructor
89  for externally allocated array, the pointer @a data can be NULL. The data
90  array can be replaced later using the method SetData().
91  */
93  : Vector(data, f->GetVSize())
94  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
95 
96  /** @brief Construct a GridFunction using previously allocated Vector @a base
97  starting at the given offset, @a base_offset. */
98  GridFunction(FiniteElementSpace *f, Vector &base, int base_offset = 0)
99  : Vector(base, base_offset, f->GetVSize())
100  { fes = f; fec = NULL; fes_sequence = f->GetSequence(); UseDevice(true); }
101 
102  /// Construct a GridFunction on the given Mesh, using the data from @a input.
103  /** The content of @a input should be in the format created by the method
104  Save(). The reconstructed FiniteElementSpace and FiniteElementCollection
105  are owned by the GridFunction. */
106  GridFunction(Mesh *m, std::istream &input);
107 
108  GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces);
109 
110  /// Copy assignment. Only the data of the base class Vector is copied.
111  /** It is assumed that this object and @a rhs use FiniteElementSpace%s that
112  have the same size.
113 
114  @note Defining this method overwrites the implicitly defined copy
115  assignment operator. */
117  { return operator=((const Vector &)rhs); }
118 
119  /// Make the GridFunction the owner of #fec and #fes.
120  /** If the new FiniteElementCollection, @a fec_, is NULL, ownership of #fec
121  and #fes is taken away. */
122  void MakeOwner(FiniteElementCollection *fec_) { fec = fec_; }
123 
125 
126  int VectorDim() const;
127  int CurlDim() const;
128 
129  /// Read only access to the (optional) internal true-dof Vector.
130  const Vector &GetTrueVector() const
131  {
132  MFEM_VERIFY(t_vec.Size() > 0, "SetTrueVector() before GetTrueVector()");
133  return t_vec;
134  }
135  /// Read and write access to the (optional) internal true-dof Vector.
136  /** Note that @a t_vec is set if it is not allocated or set already.*/
138  { if (t_vec.Size() == 0) { SetTrueVector(); } 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  /** @brief Gradient of a scalar function at a quadrature point.
331 
332  @note It is assumed that the IntegrationPoint of interest has been
333  specified by ElementTransformation::SetIntPoint() before calling
334  GetGradient().
335 
336  @note Can be used from a ParGridFunction when @a tr is an
337  ElementTransformation of a face-neighbor element and face-neighbor data
338  has been exchanged. */
339  void GetGradient(ElementTransformation &tr, Vector &grad) const;
340 
341  /// Extension of GetGradient(...) for a collection of IntegrationPoints.
343  DenseMatrix &grad) const;
344 
345  /// Extension of GetGradient(...) for a collection of IntegrationPoints.
346  void GetGradients(const int elem, const IntegrationRule &ir,
347  DenseMatrix &grad) const
348  { GetGradients(*fes->GetElementTransformation(elem), ir, grad); }
349 
350  /** @brief Compute the vector gradient with respect to the physical element
351  variable. */
352  void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const;
353 
354  /** @brief Compute the vector gradient with respect to the reference element
355  variable. */
357 
358  /** Compute \f$ (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) \f$,
359  where \f$ \psi_i \f$ are the basis functions for the FE space of avgs.
360  Both FE spaces should be scalar and on the same mesh. */
361  void GetElementAverages(GridFunction &avgs) const;
362 
363  /** Sets the output vector @a dof_vals to the values of the degrees of
364  freedom of element @a el. */
365  virtual void GetElementDofValues(int el, Vector &dof_vals) const;
366 
367  /** Impose the given bounds on the function's DOFs while preserving its local
368  * integral (described in terms of the given weights) on the i'th element
369  * through SLBPQ optimization.
370  * Intended to be used for discontinuous FE functions. */
371  void ImposeBounds(int i, const Vector &weights,
372  const Vector &lo_, const Vector &hi_);
373  void ImposeBounds(int i, const Vector &weights,
374  double min_ = 0.0, double max_ = infinity());
375 
376  /** On a non-conforming mesh, make sure the function lies in the conforming
377  space by multiplying with R and then with P, the conforming restriction
378  and prolongation matrices of the space, respectively. */
379  void RestrictConforming();
380 
381  /** @brief Project the @a src GridFunction to @a this GridFunction, both of
382  which must be on the same mesh. */
383  /** The current implementation assumes that all elements use the same
384  projection matrix. */
385  void ProjectGridFunction(const GridFunction &src);
386 
387  /** @brief Project @a coeff Coefficient to @a this GridFunction. The
388  projection computation depends on the choice of the FiniteElementSpace
389  #fes. Note that this is usually interpolation at the degrees of freedom
390  in each element (not L2 projection). */
391  virtual void ProjectCoefficient(Coefficient &coeff);
392 
393  /** @brief Project @a coeff Coefficient to @a this GridFunction, using one
394  element for each degree of freedom in @a dofs and nodal interpolation on
395  that element. */
396  void ProjectCoefficient(Coefficient &coeff, Array<int> &dofs, int vd = 0);
397 
398  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction. The
399  projection computation depends on the choice of the FiniteElementSpace
400  #fes. Note that this is usually interpolation at the degrees of freedom
401  in each element (not L2 projection).*/
403 
404  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction, using
405  one element for each degree of freedom in @a dofs and nodal interpolation
406  on that element. */
407  void ProjectCoefficient(VectorCoefficient &vcoeff, Array<int> &dofs);
408 
409  /** @brief Project @a vcoeff VectorCoefficient to @a this GridFunction, only
410  projecting onto elements with the given @a attribute */
411  void ProjectCoefficient(VectorCoefficient &vcoeff, int attribute);
412 
413  /** @brief Analogous to the version with argument @a vcoeff VectorCoefficient
414  but using an array of scalar coefficients for each component. */
415  void ProjectCoefficient(Coefficient *coeff[]);
416 
417  /** @brief Project a discontinuous vector coefficient as a grid function on
418  a continuous finite element space. The values in shared dofs are
419  determined from the element with maximal attribute. */
420  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
421 
423  /** @brief Projects a discontinuous coefficient so that the values in shared
424  vdofs are computed by taking an average of the possible values. */
425  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
426  /** @brief Projects a discontinuous _vector_ coefficient so that the values
427  in shared vdofs are computed by taking an average of the possible values.
428  */
429  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff, AvgType type);
430 
431 protected:
432  /** @brief Accumulates (depending on @a type) the values of @a coeff at all
433  shared vdofs and counts in how many zones each vdof appears. */
434  void AccumulateAndCountZones(Coefficient &coeff, AvgType type,
435  Array<int> &zones_per_vdof);
436 
437  /** @brief Accumulates (depending on @a type) the values of @a vcoeff at all
438  shared vdofs and counts in how many zones each vdof appears. */
440  Array<int> &zones_per_vdof);
441 
442  /** @brief Used for the serial and parallel implementations of the
443  GetDerivative() method; see its documentation. */
444  void AccumulateAndCountDerivativeValues(int comp, int der_comp,
445  GridFunction &der,
446  Array<int> &zones_per_dof);
447 
449  VectorCoefficient *vcoeff, Array<int> &attr,
450  Array<int> &values_counter);
451 
453  Array<int> &bdr_attr,
454  Array<int> &values_counter);
455 
456  // Complete the computation of averages; called e.g. after
457  // AccumulateAndCountZones().
458  void ComputeMeans(AvgType type, Array<int> &zones_per_vdof);
459 
460 public:
461  /** @brief For each vdof, counts how many elements contain the vdof,
462  as containment is determined by FiniteElementSpace::GetElementVDofs(). */
463  virtual void CountElementsPerVDof(Array<int> &elem_per_vdof) const;
464 
465  /** @brief Project a Coefficient on the GridFunction, modifying only DOFs on
466  the boundary associated with the boundary attributes marked in the
467  @a attr array. */
469  {
470  Coefficient *coeff_p = &coeff;
471  ProjectBdrCoefficient(&coeff_p, attr);
472  }
473 
474  /** @brief Project a VectorCoefficient on the GridFunction, modifying only
475  DOFs on the boundary associated with the boundary attributes marked in
476  the @a attr array. */
477  virtual void ProjectBdrCoefficient(VectorCoefficient &vcoeff,
478  Array<int> &attr);
479 
480  /** @brief Project a set of Coefficient%s on the components of the
481  GridFunction, modifying only DOFs on the boundary associated with the
482  boundary attributed marked in the @a attr array. */
483  /** If a Coefficient pointer in the array @a coeff is NULL, that component
484  will not be touched. */
485  virtual void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
486 
487  /** Project the normal component of the given VectorCoefficient on
488  the boundary. Only boundary attributes that are marked in
489  'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
491  Array<int> &bdr_attr);
492 
493  /** @brief Project the tangential components of the given VectorCoefficient
494  on the boundary. Only boundary attributes that are marked in @a bdr_attr
495  are projected. Assumes ND-type VectorFE GridFunction. */
496  virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff,
497  Array<int> &bdr_attr);
498 
499 
500  virtual double ComputeL2Error(Coefficient *exsol[],
501  const IntegrationRule *irs[] = NULL,
502  const Array<int> *elems = NULL) const;
503 
504  /// Returns ||grad u_ex - grad u_h||_L2 in element ielem for H1 or L2 elements
505  virtual double ComputeElementGradError(int ielem, VectorCoefficient *exgrad,
506  const IntegrationRule *irs[] = NULL) const;
507 
508  /// Returns ||u_ex - u_h||_L2 for H1 or L2 elements
509  /* The @a elems input variable expects a list of markers:
510  an elem marker equal to 1 will compute the L2 error on that element
511  an elem marker equal to 0 will not compute the L2 error on that element */
512  virtual double ComputeL2Error(Coefficient &exsol,
513  const IntegrationRule *irs[] = NULL,
514  const Array<int> *elems = NULL) const
515  { return GridFunction::ComputeLpError(2.0, exsol, NULL, irs, elems); }
516 
517  virtual double ComputeL2Error(VectorCoefficient &exsol,
518  const IntegrationRule *irs[] = NULL,
519  const Array<int> *elems = NULL) const;
520 
521  /// Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements
522  virtual double ComputeGradError(VectorCoefficient *exgrad,
523  const IntegrationRule *irs[] = NULL) const;
524 
525  /// Returns ||curl u_ex - curl u_h||_L2 for ND elements
526  virtual double ComputeCurlError(VectorCoefficient *excurl,
527  const IntegrationRule *irs[] = NULL) const;
528 
529  /// Returns ||div u_ex - div u_h||_L2 for RT elements
530  virtual double ComputeDivError(Coefficient *exdiv,
531  const IntegrationRule *irs[] = NULL) const;
532 
533  /// Returns the Face Jumps error for L2 elements. The error can be weighted
534  /// by a constant nu, by nu/h, or nu*p^2/h, depending on the value of
535  /// @a jump_scaling.
536  virtual double ComputeDGFaceJumpError(Coefficient *exsol,
537  Coefficient *ell_coeff,
538  class JumpScaling jump_scaling,
539  const IntegrationRule *irs[] = NULL)
540  const;
541 
542  /// Returns the Face Jumps error for L2 elements, with 1/h scaling.
543  MFEM_DEPRECATED
544  double ComputeDGFaceJumpError(Coefficient *exsol,
545  Coefficient *ell_coeff,
546  double Nu,
547  const IntegrationRule *irs[] = NULL) const;
548 
549  /** This method is kept for backward compatibility.
550 
551  Returns either the H1-seminorm, or the DG face jumps error, or both
552  depending on norm_type = 1, 2, 3. Additional arguments for the DG face
553  jumps norm: ell_coeff: mesh-depended coefficient (weight) Nu: scalar
554  constant weight */
555  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
556  Coefficient *ell_coef, double Nu,
557  int norm_type) const;
558 
559  /// Returns the error measured in H1-norm for H1 elements or in "broken"
560  /// H1-norm for L2 elements
561  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
562  const IntegrationRule *irs[] = NULL) const;
563 
564  /// Returns the error measured in H(div)-norm for RT elements
565  virtual double ComputeHDivError(VectorCoefficient *exsol,
566  Coefficient *exdiv,
567  const IntegrationRule *irs[] = NULL) const;
568 
569  /// Returns the error measured in H(curl)-norm for ND elements
570  virtual double ComputeHCurlError(VectorCoefficient *exsol,
571  VectorCoefficient *excurl,
572  const IntegrationRule *irs[] = NULL) const;
573 
574  virtual double ComputeMaxError(Coefficient &exsol,
575  const IntegrationRule *irs[] = NULL) const
576  {
577  return ComputeLpError(infinity(), exsol, NULL, irs);
578  }
579 
580  virtual double ComputeMaxError(Coefficient *exsol[],
581  const IntegrationRule *irs[] = NULL) const;
582 
583  virtual double ComputeMaxError(VectorCoefficient &exsol,
584  const IntegrationRule *irs[] = NULL) const
585  {
586  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
587  }
588 
589  virtual double ComputeL1Error(Coefficient &exsol,
590  const IntegrationRule *irs[] = NULL) const
591  { return ComputeLpError(1.0, exsol, NULL, irs); }
592 
593  virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
594  int norm_type, const Array<int> *elems = NULL,
595  const IntegrationRule *irs[] = NULL) const;
596 
597  virtual double ComputeL1Error(VectorCoefficient &exsol,
598  const IntegrationRule *irs[] = NULL) const
599  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
600 
601  /* The @a elems input variable expects a list of markers:
602  an elem marker equal to 1 will compute the L2 error on that element
603  an elem marker equal to 0 will not compute the L2 error on that element */
604  virtual double ComputeLpError(const double p, Coefficient &exsol,
605  Coefficient *weight = NULL,
606  const IntegrationRule *irs[] = NULL,
607  const Array<int> *elems = NULL) const;
608 
609  /** Compute the Lp error in each element of the mesh and store the results in
610  the Vector @a 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, Coefficient &exsol,
613  Vector &error,
614  Coefficient *weight = NULL,
615  const IntegrationRule *irs[] = NULL
616  ) const;
617 
618  virtual void ComputeElementL1Errors(Coefficient &exsol,
619  Vector &error,
620  const IntegrationRule *irs[] = NULL
621  ) const
622  { ComputeElementLpErrors(1.0, exsol, error, NULL, irs); }
623 
624  virtual void ComputeElementL2Errors(Coefficient &exsol,
625  Vector &error,
626  const IntegrationRule *irs[] = NULL
627  ) const
628  { ComputeElementLpErrors(2.0, exsol, error, NULL, irs); }
629 
631  Vector &error,
632  const IntegrationRule *irs[] = NULL
633  ) const
634  { ComputeElementLpErrors(infinity(), exsol, error, NULL, irs); }
635 
636  /** When given a vector weight, compute the pointwise (scalar) error as the
637  dot product of the vector error with the vector weight. Otherwise, the
638  scalar error is the l_2 norm of the vector error. */
639  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
640  Coefficient *weight = NULL,
641  VectorCoefficient *v_weight = NULL,
642  const IntegrationRule *irs[] = NULL) const;
643 
644  /** Compute the Lp error in each element of the mesh and store the results in
645  the Vector @ error. The result should be of length number of elements,
646  for example an L2 GridFunction of order zero using map type VALUE. */
647  virtual void ComputeElementLpErrors(const double p, VectorCoefficient &exsol,
648  Vector &error,
649  Coefficient *weight = NULL,
650  VectorCoefficient *v_weight = NULL,
651  const IntegrationRule *irs[] = NULL
652  ) const;
653 
655  Vector &error,
656  const IntegrationRule *irs[] = NULL
657  ) const
658  { ComputeElementLpErrors(1.0, exsol, error, NULL, NULL, irs); }
659 
661  Vector &error,
662  const IntegrationRule *irs[] = NULL
663  ) const
664  { ComputeElementLpErrors(2.0, exsol, error, NULL, NULL, irs); }
665 
667  Vector &error,
668  const IntegrationRule *irs[] = NULL
669  ) const
670  { ComputeElementLpErrors(infinity(), exsol, error, NULL, NULL, irs); }
671 
672  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
673  GridFunction &flux,
674  bool wcoef = true, int subdomain = -1);
675 
676  /// Redefine '=' for GridFunction = constant.
677  GridFunction &operator=(double value);
678 
679  /// Copy the data from @a v.
680  /** The size of @a v must be equal to the size of the associated
681  FiniteElementSpace #fes. */
682  GridFunction &operator=(const Vector &v);
683 
684  /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
685  virtual void Update();
686 
687  /** Return update counter, similar to Mesh::GetSequence(). Used to
688  check if it is up to date with the space. */
689  long GetSequence() const { return fes_sequence; }
690 
692  const FiniteElementSpace *FESpace() const { return fes; }
693 
694  /// Associate a new FiniteElementSpace with the GridFunction.
695  /** The GridFunction is resized using the SetSize() method. */
696  virtual void SetSpace(FiniteElementSpace *f);
697 
698  using Vector::MakeRef;
699 
700  /** @brief Make the GridFunction reference external data on a new
701  FiniteElementSpace. */
702  /** This method changes the FiniteElementSpace associated with the
703  GridFunction and sets the pointer @a v as external data in the
704  GridFunction. */
705  virtual void MakeRef(FiniteElementSpace *f, double *v);
706 
707  /** @brief Make the GridFunction reference external data on a new
708  FiniteElementSpace. */
709  /** This method changes the FiniteElementSpace associated with the
710  GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
711  as external data in the GridFunction.
712  @note This version of the method will also perform bounds checks when
713  the build option MFEM_DEBUG is enabled. */
714  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
715 
716  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
717  GridFunction. */
718  /** - If the prolongation matrix of @a f is trivial (i.e. its method
719  FiniteElementSpace::GetProlongationMatrix() returns NULL), then the
720  method MakeRef() is called with the same arguments.
721  - Otherwise, the method SetSpace() is called with argument @a f.
722  - The internal true-dof vector is set to reference @a tv. */
723  void MakeTRef(FiniteElementSpace *f, double *tv);
724 
725  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
726  GridFunction. */
727  /** - If the prolongation matrix of @a f is trivial (i.e. its method
728  FiniteElementSpace::GetProlongationMatrix() returns NULL), this method
729  calls MakeRef() with the same arguments.
730  - Otherwise, this method calls SetSpace() with argument @a f.
731  - The internal true-dof vector is set to reference the sub-vector of
732  @a tv starting at the offset @a tv_offset. */
733  void MakeTRef(FiniteElementSpace *f, Vector &tv, int tv_offset);
734 
735  /// Save the GridFunction to an output stream.
736  virtual void Save(std::ostream &out) const;
737 
738  /// Save the GridFunction to a file. The given @a precision will be used for
739  /// ASCII output.
740  virtual void Save(const char *fname, int precision=16) const;
741 
742 #ifdef MFEM_USE_ADIOS2
743  /// Save the GridFunction to a binary output stream using adios2 bp format.
744  virtual void Save(adios2stream &out, const std::string& variable_name,
747 #endif
748 
749  /** @brief Write the GridFunction in VTK format. Note that Mesh::PrintVTK
750  must be called first. The parameter ref > 0 must match the one used in
751  Mesh::PrintVTK. */
752  void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
753 
754  /** @brief Write the GridFunction in STL format. Note that the mesh dimension
755  must be 2 and that quad elements will be broken into two triangles.*/
756  void SaveSTL(std::ostream &out, int TimesToRefine = 1);
757 
758  /// Destroys grid function.
759  virtual ~GridFunction() { Destroy(); }
760 };
761 
762 
763 /** Overload operator<< for std::ostream and GridFunction; valid also for the
764  derived class ParGridFunction */
765 std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
766 
767 /// Class used to specify how the jump terms in
768 /// GridFunction::ComputeDGFaceJumpError are scaled.
770 {
771 public:
773  {
777  };
778 private:
779  double nu;
780  JumpScalingType type;
781 public:
782  JumpScaling(double nu_=1.0, JumpScalingType type_=CONSTANT)
783  : nu(nu_), type(type_) { }
784  double Eval(double h, int p) const
785  {
786  double val = nu;
787  if (type != CONSTANT) { val /= h; }
788  if (type == P_SQUARED_OVER_H) { val *= p*p; }
789  return val;
790  }
791 };
792 
793 /// Overload operator<< for std::ostream and QuadratureFunction.
794 std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
795 
796 
797 double ZZErrorEstimator(BilinearFormIntegrator &blfi,
798  GridFunction &u,
799  GridFunction &flux,
800  Vector &error_estimates,
801  Array<int> *aniso_flags = NULL,
802  int with_subdomains = 1,
803  bool with_coeff = false);
804 
805 /// Defines the global tensor product polynomial space used by NewZZErorrEstimator
806 /**
807  * See BoundingBox(...) for a description of @a angle and @a midpoint
808  */
809 void TensorProductLegendre(int dim, // input
810  int order, // input
811  const Vector &x_in, // input
812  const Vector &xmax, // input
813  const Vector &xmin, // input
814  Vector &poly, // output
815  double angle=0.0, // input (optional)
816  const Vector *midpoint=NULL); // input (optional)
817 
818 /// Defines the bounding box for the face patches used by NewZZErorrEstimator
819 /**
820  * By default, BoundingBox(...) computes the parameters of a minimal bounding box
821  * for the given @a face_patch that is aligned with the physical (i.e. global)
822  * Cartesian axes. This means that the size of the bounding box will depend on the
823  * orientation of the patch. It is better to construct an orientation-independent box.
824  * This is implemented for 2D patches. The parameters @a angle and @a midpoint encode
825  * the necessary additional geometric information.
826  *
827  * @a iface : Index of the face that the patch corresponds to.
828  * This is used to compute @a angle and @a midpoint.
829  *
830  * @a angle : The angle the patch face makes with the x-axis.
831  * @a midpoint : The midpoint of the face.
832  */
833 void BoundingBox(const Array<int> &face_patch, // input
834  FiniteElementSpace *ufes, // input
835  int order, // input
836  Vector &xmin, // output
837  Vector &xmax, // output
838  double &angle, // output
839  Vector &midpoint, // output
840  int iface=-1); // input (optional)
841 
842 /// A ``true'' ZZ error estimator that uses face-based patches for flux reconstruction.
843 /**
844  * Only two-element face patches are ever used:
845  * - For conforming faces, the face patch consists of its two neighboring elements.
846  * - In the non-conforming setting, only the face patches associated to fine-scale
847  * element faces are used. These face patches always consist of two elements
848  * delivered by mesh::GetFaceElements(Face, *Elem1, *Elem2).
849  */
850 double LSZZErrorEstimator(BilinearFormIntegrator &blfi, // input
851  GridFunction &u, // input
852  Vector &error_estimates, // output
853  bool subdomain_reconstruction = true, // input (optional)
854  bool with_coeff = false, // input (optional)
855  double tichonov_coeff = 0.0); // input (optional)
856 
857 /// Compute the Lp distance between two grid functions on the given element.
858 double ComputeElementLpDistance(double p, int i,
859  GridFunction& gf1, GridFunction& gf2);
860 
861 
862 /// Class used for extruding scalar GridFunctions
864 {
865 private:
866  int n;
867  Mesh *mesh_in;
868  Coefficient &sol_in;
869 public:
871  : n(n_), mesh_in(m), sol_in(s) { }
872  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
873  virtual ~ExtrudeCoefficient() { }
874 };
875 
876 /// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
877 GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
878  GridFunction *sol, const int ny);
879 
880 } // namespace mfem
881 
882 #endif
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:3816
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:2786
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.hpp:759
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:2920
virtual void CountElementsPerVDof(Array< int > &elem_per_vdof) const
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
Definition: gridfunc.cpp:2016
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:96
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Memory< double > data
Definition: vector.hpp:62
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:863
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:476
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
Definition: gridfunc.cpp:251
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.
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:2348
void RestrictConforming()
Definition: gridfunc.cpp:1976
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:4529
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
Definition: gridfunc.cpp:1736
int CurlDim() const
Definition: gridfunc.cpp:346
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:2589
virtual void ComputeElementL1Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:654
virtual void ComputeElementLpErrors(const double p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3433
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
Definition: gridfunc.cpp:1643
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
void GetElementAverages(GridFunction &avgs) const
Definition: gridfunc.cpp:1798
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 ...
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, Array< int > &bdr_attr, Array< int > &values_counter)
Definition: gridfunc.cpp:2271
long GetSequence() const
Definition: gridfunc.hpp:689
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:137
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:122
double Eval(double h, int p) const
Definition: gridfunc.hpp:784
ExtrudeCoefficient(Mesh *m, Coefficient &s, int n_)
Definition: gridfunc.hpp:870
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:2323
virtual void ComputeElementMaxErrors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:630
void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess, int vdim=1) const
Definition: gridfunc.cpp:608
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:449
void ImposeBounds(int i, const Vector &weights, const Vector &lo_, const Vector &hi_)
Definition: gridfunc.cpp:1901
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:597
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Definition: gridfunc.cpp:4519
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:574
GridFunction(FiniteElementSpace *f, double *data)
Construct a GridFunction using previously allocated array data.
Definition: gridfunc.hpp:92
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:1430
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition: gridfunc.cpp:366
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:2694
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:144
void GetCurl(ElementTransformation &tr, Vector &curl) const
Definition: gridfunc.cpp:1547
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:521
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
Definition: gridfunc.hpp:346
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:583
GridFunction(const GridFunction &orig)
Copy constructor. The internal true-dof vector t_vec is not copied.
Definition: gridfunc.hpp:78
JumpScaling(double nu_=1.0, JumpScalingType type_=CONSTANT)
Definition: gridfunc.hpp:782
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:221
void BoundingBox(const Array< int > &patch, FiniteElementSpace *ufes, int order, Vector &xmin, Vector &xmax, double &angle, Vector &midpoint, int iface)
Defines the bounding box for the face patches used by NewZZErorrEstimator.
Definition: gridfunc.cpp:4188
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:3186
virtual double ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
Definition: gridfunc.cpp:2960
double LSZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, Vector &error_estimates, bool subdomain_reconstruction, bool with_coeff, double tichonov_coeff)
A ‘‘true’’ ZZ error estimator that uses face-based patches for flux reconstruction.
Definition: gridfunc.cpp:4255
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:1377
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:2033
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:203
double GetDivergence(ElementTransformation &tr) const
Definition: gridfunc.cpp:1457
FiniteElementCollection * OwnFEC()
Definition: gridfunc.hpp:124
virtual ~ExtrudeCoefficient()
Definition: gridfunc.hpp:873
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:691
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition: gridfunc.hpp:116
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3260
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:395
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:589
virtual void GetElementDofValues(int el, Vector &dof_vals) const
Definition: gridfunc.cpp:1838
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
virtual double ComputeElementGradError(int ielem, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 in element ielem for H1 or L2 elements.
Definition: gridfunc.cpp:2882
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
virtual double ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:3037
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
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
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:3160
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
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
Definition: gridfunc.cpp:1849
virtual double ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
Definition: gridfunc.cpp:3001
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:1306
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||u_ex - u_h||_L2 for H1 or L2 elements.
Definition: gridfunc.hpp:512
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 GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Compute the vector gradient with respect to the reference element variable.
Definition: gridfunc.cpp:1441
void TensorProductLegendre(int dim, int order, const Vector &x_in, const Vector &xmax, const Vector &xmin, Vector &poly, double angle, const Vector *midpoint)
Defines the global tensor product polynomial space used by NewZZErorrEstimator.
Definition: gridfunc.cpp:4096
virtual void ComputeElementL2Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:624
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:4462
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:130
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:3195
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:3366
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains, bool with_coeff)
Definition: gridfunc.cpp:3998
Class for integration point with weight.
Definition: intrules.hpp:31
int VectorDim() const
Definition: gridfunc.cpp:323
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:3836
virtual void ComputeElementL1Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:618
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:1241
int dim
Definition: ex24.cpp:53
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
Definition: gridfunc.cpp:1712
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:2769
void LegacyNCReorder()
Loading helper.
Definition: gridfunc.cpp:3936
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:664
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition: vector.hpp:44
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:1114
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:98
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2415
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:195
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:1279
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:711
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:714
void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps, int vdim=1) const
Definition: gridfunc.cpp:569
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:118
Vector data type.
Definition: vector.hpp:58
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:581
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition: gridfunc.cpp:308
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Definition: gridfunc.hpp:83
RefCoord s[3]
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3696
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: gridfunc.cpp:381
virtual void ComputeElementMaxErrors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:666
const FiniteElementSpace * FESpace() const
Definition: gridfunc.hpp:692
void GetValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:1157
void GetBdrValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:1203
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr, Array< int > &values_counter)
Definition: gridfunc.cpp:2129
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 > 0 must match the one used in Mesh::PrintVTK.
Definition: gridfunc.cpp:3737
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:468
virtual void ComputeElementL2Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:660
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:1337
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:769