MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
gridfunc.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_GRIDFUNC
13 #define MFEM_GRIDFUNC
14 
15 #include "../config/config.hpp"
16 #include "fespace.hpp"
17 #include "coefficient.hpp"
18 #include "bilininteg.hpp"
19 #include <limits>
20 #include <ostream>
21 #include <string>
22 
23 namespace mfem
24 {
25 
26 /// Class for grid function - Vector with associated FE space.
27 class GridFunction : public Vector
28 {
29 protected:
30  /// FE space on which grid function lives.
32 
33  /// Used when the grid function is read from a file
35 
36  long sequence; // see FiniteElementSpace::sequence, Mesh::sequence
37 
38  /** Optional, internal true-dof vector: if the FiniteElementSpace #fes has a
39  non-trivial (i.e. not NULL) prolongation operator, this Vector may hold
40  associated true-dof values - either owned or external. */
42 
43  void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[]);
44 
46 
47  // Project the delta coefficient without scaling and return the (local)
48  // integral of the projection.
49  void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff,
50  double &integral);
51 
52  // Sum fluxes to vertices and count element contributions
54  GridFunction &flux,
55  Array<int>& counts,
56  int wcoef,
57  int subdomain);
58 
59  /** Project a discontinuous vector coefficient in a continuous space and
60  return in dof_attr the maximal attribute of the elements containing each
61  degree of freedom. */
62  void ProjectDiscCoefficient(VectorCoefficient &coeff, Array<int> &dof_attr);
63 
64  void Destroy();
65 
66 public:
67 
68  GridFunction() { fes = NULL; fec = NULL; sequence = 0; }
69 
70  /// Copy constructor.
72  : Vector(orig), fes(orig.fes), fec(NULL), sequence(orig.sequence) { }
73 
74  /// Construct a GridFunction associated with the FiniteElementSpace @a *f.
75  GridFunction(FiniteElementSpace *f) : Vector(f->GetVSize())
76  { fes = f; fec = NULL; sequence = f->GetSequence(); }
77 
78  /// Construct a GridFunction using previously allocated array @a data.
79  /** The GridFunction does not assume ownership of @a data which is assumed to
80  be of size at least `f->GetVSize()`. Similar to the Vector constructor
81  for externally allocated array, the pointer @a data can be NULL. The data
82  array can be replaced later using the method SetData().
83  */
84  GridFunction(FiniteElementSpace *f, double *data) : Vector(data, f->GetVSize())
85  { fes = f; fec = NULL; sequence = f->GetSequence(); }
86 
87  /// Construct a GridFunction on the given Mesh, using the data from @a input.
88  /** The content of @a input should be in the format created by the method
89  Save(). The reconstructed FiniteElementSpace and FiniteElementCollection
90  are owned by the GridFunction. */
91  GridFunction(Mesh *m, std::istream &input);
92 
93  GridFunction(Mesh *m, GridFunction *gf_array[], int num_pieces);
94 
95  /// Make the GridFunction the owner of 'fec' and 'fes'
96  void MakeOwner(FiniteElementCollection *_fec) { fec = _fec; }
97 
99 
100  int VectorDim() const;
101 
102  /// Read only access to the (optional) internal true-dof Vector.
103  /** Note that the returned Vector may be empty, if not previously allocated
104  or set. */
105  const Vector &GetTrueVector() const { return t_vec; }
106  /// Read and write access to the (optional) internal true-dof Vector.
107  /** Note that the returned Vector may be empty, if not previously allocated
108  or set. */
109  Vector &GetTrueVector() { return t_vec; }
110 
111  /// @brief Extract the true-dofs from the GridFunction. If all dofs are true,
112  /// then `tv` will be set to point to the data of `*this`.
113  void GetTrueDofs(Vector &tv) const;
114 
115  /// Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
117 
118  /// Set the GridFunction from the given true-dof vector.
119  virtual void SetFromTrueDofs(const Vector &tv);
120 
121  /// Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
123 
124  /// Returns the values in the vertices of i'th element for dimension vdim.
125  void GetNodalValues(int i, Array<double> &nval, int vdim = 1) const;
126 
127  virtual double GetValue(int i, const IntegrationPoint &ip,
128  int vdim = 1) const;
129 
130  void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const;
131 
132  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
133  int vdim = 1) const;
134 
135  void GetValues(int i, const IntegrationRule &ir, Vector &vals,
136  DenseMatrix &tr, int vdim = 1) const;
137 
138  int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals,
139  DenseMatrix &tr, int vdim = 1) const;
140 
142  DenseMatrix &vals) const;
143 
144  void GetVectorValues(int i, const IntegrationRule &ir,
145  DenseMatrix &vals, DenseMatrix &tr) const;
146 
147  int GetFaceVectorValues(int i, int side, const IntegrationRule &ir,
148  DenseMatrix &vals, DenseMatrix &tr) const;
149 
150  void GetValuesFrom(const GridFunction &orig_func);
151 
152  void GetBdrValuesFrom(const GridFunction &orig_func);
153 
154  void GetVectorFieldValues(int i, const IntegrationRule &ir,
155  DenseMatrix &vals,
156  DenseMatrix &tr, int comp = 0) const;
157 
158  /// For a vector grid function, makes sure that the ordering is byNODES.
159  void ReorderByNodes();
160 
161  /// Return the values as a vector on mesh vertices for dimension vdim.
162  void GetNodalValues(Vector &nval, int vdim = 1) const;
163 
164  void GetVectorFieldNodalValues(Vector &val, int comp) const;
165 
166  void ProjectVectorFieldOn(GridFunction &vec_field, int comp = 0);
167 
168  void GetDerivative(int comp, int der_comp, GridFunction &der);
169 
170  double GetDivergence(ElementTransformation &tr) const;
171 
172  void GetCurl(ElementTransformation &tr, Vector &curl) const;
173 
174  void GetGradient(ElementTransformation &tr, Vector &grad) const;
175 
176  void GetGradients(const int elem, const IntegrationRule &ir,
177  DenseMatrix &grad) const;
178 
179  void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const;
180 
181  /** Compute \f$ (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) \f$,
182  where \f$ \psi_i \f$ are the basis functions for the FE space of avgs.
183  Both FE spaces should be scalar and on the same mesh. */
184  void GetElementAverages(GridFunction &avgs) const;
185 
186  /** Impose the given bounds on the function's DOFs while preserving its local
187  * integral (described in terms of the given weights) on the i'th element
188  * through SLBPQ optimization.
189  * Intended to be used for discontinuous FE functions. */
190  void ImposeBounds(int i, const Vector &weights,
191  const Vector &_lo, const Vector &_hi);
192  void ImposeBounds(int i, const Vector &weights,
193  double _min = 0.0, double _max = infinity());
194 
195  /** @brief Project the @a src GridFunction to @a this GridFunction, both of
196  which must be on the same mesh. */
197  /** The current implementation assumes that all elements use the same
198  projection matrix. */
199  void ProjectGridFunction(const GridFunction &src);
200 
201  virtual void ProjectCoefficient(Coefficient &coeff);
202 
203  // call fes -> BuildDofToArrays() before using this projection
204  void ProjectCoefficient(Coefficient &coeff, Array<int> &dofs, int vd = 0);
205 
207 
208  // call fes -> BuildDofToArrays() before using this projection
209  void ProjectCoefficient(VectorCoefficient &vcoeff, Array<int> &dofs);
210 
211  void ProjectCoefficient(Coefficient *coeff[]);
212 
213  /** @brief Project a discontinuous vector coefficient as a grid function on
214  a continuous finite element space. The values in shared dofs are
215  determined from the element with maximal attribute. */
216  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff);
217 
219  /** @brief Projects a discontinuous coefficient so that the values in shared
220  vdofs are computed by taking an average of the possible values. */
221  virtual void ProjectDiscCoefficient(Coefficient &coeff, AvgType type);
222  /** @brief Projects a discontinuous _vector_ coefficient so that the values
223  in shared vdofs are computed by taking an average of the possible values.
224  */
225  virtual void ProjectDiscCoefficient(VectorCoefficient &coeff, AvgType type);
226 
227 protected:
228  /** @brief Accumulates (depending on @a type) the values of @a coeff at all
229  shared vdofs and counts in how many zones each vdof appears. */
230  void AccumulateAndCountZones(Coefficient &coeff, AvgType type,
231  Array<int> &zones_per_vdof);
232 
233  /** @brief Accumulates (depending on @a type) the values of @a vcoeff at all
234  shared vdofs and counts in how many zones each vdof appears. */
236  Array<int> &zones_per_vdof);
237 
238  // Complete the computation of averages; called e.g. after
239  // AccumulateAndCountZones().
240  void ComputeMeans(AvgType type, Array<int> &zones_per_vdof);
241 
242 public:
244  {
245  Coefficient *coeff_p = &coeff;
246  ProjectBdrCoefficient(&coeff_p, attr);
247  }
248 
249  void ProjectBdrCoefficient(Coefficient *coeff[], Array<int> &attr);
250 
251  /** Project the normal component of the given VectorCoefficient on
252  the boundary. Only boundary attributes that are marked in
253  'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
255  Array<int> &bdr_attr);
256 
257  /** Project the tangential components of the given VectorCoefficient on
258  the boundary. Only boundary attributes that are marked in
259  'bdr_attr' are projected. Assumes ND-type VectorFE GridFunction. */
261  Array<int> &bdr_attr);
262 
263  virtual double ComputeL2Error(Coefficient &exsol,
264  const IntegrationRule *irs[] = NULL) const
265  { return ComputeLpError(2.0, exsol, NULL, irs); }
266 
267  virtual double ComputeL2Error(Coefficient *exsol[],
268  const IntegrationRule *irs[] = NULL) const;
269 
270  virtual double ComputeL2Error(VectorCoefficient &exsol,
271  const IntegrationRule *irs[] = NULL,
272  Array<int> *elems = NULL) const;
273 
274  virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
275  Coefficient *ell_coef, double Nu,
276  int norm_type) const;
277 
278  virtual double ComputeMaxError(Coefficient &exsol,
279  const IntegrationRule *irs[] = NULL) const
280  {
281  return ComputeLpError(infinity(), exsol, NULL, irs);
282  }
283 
284  virtual double ComputeMaxError(Coefficient *exsol[],
285  const IntegrationRule *irs[] = NULL) const;
286 
287  virtual double ComputeMaxError(VectorCoefficient &exsol,
288  const IntegrationRule *irs[] = NULL) const
289  {
290  return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
291  }
292 
293  virtual double ComputeL1Error(Coefficient &exsol,
294  const IntegrationRule *irs[] = NULL) const
295  { return ComputeLpError(1.0, exsol, NULL, irs); }
296 
297  virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
298  int norm_type, Array<int> *elems = NULL,
299  const IntegrationRule *irs[] = NULL) const;
300 
301  virtual double ComputeL1Error(VectorCoefficient &exsol,
302  const IntegrationRule *irs[] = NULL) const
303  { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
304 
305  virtual double ComputeLpError(const double p, Coefficient &exsol,
306  Coefficient *weight = NULL,
307  const IntegrationRule *irs[] = NULL) const;
308 
309  /** When given a vector weight, compute the pointwise (scalar) error as the
310  dot product of the vector error with the vector weight. Otherwise, the
311  scalar error is the l_2 norm of the vector error. */
312  virtual double ComputeLpError(const double p, VectorCoefficient &exsol,
313  Coefficient *weight = NULL,
314  VectorCoefficient *v_weight = NULL,
315  const IntegrationRule *irs[] = NULL) const;
316 
317  virtual void ComputeFlux(BilinearFormIntegrator &blfi,
318  GridFunction &flux,
319  int wcoef = 1, int subdomain = -1);
320 
321  /// Redefine '=' for GridFunction = constant.
322  GridFunction &operator=(double value);
323 
324  /// Copy the data from @a v.
325  /** The size of @a v must be equal to the size of the FiniteElementSpace
326  @a fes. */
327  GridFunction &operator=(const Vector &v);
328 
329  /// Copy the data from @a v.
330  /** The GridFunctions @a v and @a *this must have FiniteElementSpaces with
331  the same size. */
333 
334  /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
335  virtual void Update();
336 
338  const FiniteElementSpace *FESpace() const { return fes; }
339 
340  /// Associate a new FiniteElementSpace with the GridFunction.
341  /** The GridFunction is resized using the SetSize() method. */
342  virtual void SetSpace(FiniteElementSpace *f);
343 
344  /** @brief Make the GridFunction reference external data on a new
345  FiniteElementSpace. */
346  /** This method changes the FiniteElementSpace associated with the
347  GridFunction and sets the pointer @a v as external data in the
348  GridFunction. */
349  virtual void MakeRef(FiniteElementSpace *f, double *v);
350 
351  /** @brief Make the GridFunction reference external data on a new
352  FiniteElementSpace. */
353  /** This method changes the FiniteElementSpace associated with the
354  GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
355  as external data in the GridFunction.
356  @note This version of the method will also perform bounds checks when
357  the build option MFEM_DEBUG is enabled. */
358  virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
359 
360  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
361  GridFunction. */
362  /** - If the prolongation matrix of @a f is trivial (i.e. its method
363  FiniteElementSpace::GetProlongationMatrix() returns NULL), then the
364  method MakeRef() is called with the same arguments.
365  - Otherwise, the method SetSpace() is called with argument @a f.
366  - The internal true-dof vector is set to reference @a tv. */
367  void MakeTRef(FiniteElementSpace *f, double *tv);
368 
369  /** @brief Associate a new FiniteElementSpace and new true-dof data with the
370  GridFunction. */
371  /** - If the prolongation matrix of @a f is trivial (i.e. its method
372  FiniteElementSpace::GetProlongationMatrix() returns NULL), this method
373  calls MakeRef() with the same arguments.
374  - Otherwise, this method calls SetSpace() with argument @a f.
375  - The internal true-dof vector is set to reference the sub-vector of
376  @a tv starting at the offset @a tv_offset. */
377  void MakeTRef(FiniteElementSpace *f, Vector &tv, int tv_offset);
378 
379  /// Save the GridFunction to an output stream.
380  virtual void Save(std::ostream &out) const;
381 
382  /** Write the GridFunction in VTK format. Note that Mesh::PrintVTK must be
383  called first. The parameter ref > 0 must match the one used in
384  Mesh::PrintVTK. */
385  void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
386 
387  void SaveSTL(std::ostream &out, int TimesToRefine = 1);
388 
389  /// Destroys grid function.
390  virtual ~GridFunction() { Destroy(); }
391 };
392 
393 
394 /** Overload operator<< for std::ostream and GridFunction; valid also for the
395  derived class ParGridFunction */
396 std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
397 
398 
399 /** @brief Class representing a function through its values (scalar or vector)
400  at quadrature points. */
402 {
403 protected:
404  QuadratureSpace *qspace; ///< Associated QuadratureSpace
405  int vdim; ///< Vector dimension
406  bool own_qspace; ///< QuadratureSpace ownership flag
407 
408 public:
409  /// Create an empty QuadratureFunction.
410  /** The object can be initialized later using the SetSpace() methods. */
412  : qspace(NULL), vdim(0), own_qspace(false) { }
413 
414  /// Create a QuadratureFunction based on the given QuadratureSpace.
415  /** The QuadratureFunction does not assume ownership of the QuadratureSpace.
416  @note The Vector data is not initialized. */
417  QuadratureFunction(QuadratureSpace *qspace_, int vdim_ = 1)
418  : Vector(vdim_*qspace_->GetSize()),
419  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
420 
421  /** @brief Create a QuadratureFunction based on the given QuadratureSpace,
422  using the external data, @a qf_data. */
423  /** The QuadratureFunction does not assume ownership of neither the
424  QuadratureSpace nor the external data. */
425  QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_ = 1)
426  : Vector(qf_data, vdim_*qspace_->GetSize()),
427  qspace(qspace_), vdim(vdim_), own_qspace(false) { }
428 
429  /// Read a QuadratureFunction from the stream @a in.
430  /** The QuadratureFunction assumes ownership of the read QuadratureSpace. */
431  QuadratureFunction(Mesh *mesh, std::istream &in);
432 
433  virtual ~QuadratureFunction() { if (own_qspace) { delete qspace; } }
434 
435  /// Get the associated QuadratureSpace.
436  QuadratureSpace *GetSpace() const { return qspace; }
437 
438  /// Change the QuadratureSpace and optionally the vector dimension.
439  /** If the new QuadratureSpace is different from the current one, the
440  QuadratureFunction will not assume ownership of the new space; otherwise,
441  the ownership flag remains the same.
442 
443  If the new vector dimension @a vdim_ < 0, the vector dimension remains
444  the same.
445 
446  The data size is updated by calling Vector::SetSize(). */
447  inline void SetSpace(QuadratureSpace *qspace_, int vdim_ = -1);
448 
449  /** @brief Change the QuadratureSpace, the data array, and optionally the
450  vector dimension. */
451  /** If the new QuadratureSpace is different from the current one, the
452  QuadratureFunction will not assume ownership of the new space; otherwise,
453  the ownership flag remains the same.
454 
455  If the new vector dimension @a vdim_ < 0, the vector dimension remains
456  the same.
457 
458  The data array is replaced by calling Vector::NewDataAndSize(). */
459  inline void SetSpace(QuadratureSpace *qspace_, double *qf_data,
460  int vdim_ = -1);
461 
462  /// Get the vector dimension.
463  int GetVDim() const { return vdim; }
464 
465  /// Set the vector dimension, updating the size by calling Vector::SetSize().
466  void SetVDim(int vdim_)
467  { vdim = vdim_; SetSize(vdim*qspace->GetSize()); }
468 
469  /// Get the QuadratureSpace ownership flag.
470  bool OwnsSpace() { return own_qspace; }
471 
472  /// Set the QuadratureSpace ownership flag.
473  void SetOwnsSpace(bool own) { own_qspace = own; }
474 
475  /// Redefine '=' for QuadratureFunction = constant.
476  QuadratureFunction &operator=(double value);
477 
478  /// Copy the data from @a v.
479  /** The size of @a v must be equal to the size of the QuadratureSpace
480  @a qspace. */
482 
483  /// Copy the data from @a v.
484  /** The QuadratureFunctions @a v and @a *this must have QuadratureSpaces with
485  the same size. */
487 
488  /// Get the IntegrationRule associated with mesh element @a idx.
489  const IntegrationRule &GetElementIntRule(int idx) const
490  { return qspace->GetElementIntRule(idx); }
491 
492  /// Return all values associated with mesh element @a idx in a Vector.
493  /** The result is stored in the Vector @a values as a reference to the
494  global values.
495 
496  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
497  `i`-th vector component at the `j`-th quadrature point.
498  */
499  inline void GetElementValues(int idx, Vector &values);
500 
501  /// Return all values associated with mesh element @a idx in a Vector.
502  /** The result is stored in the Vector @a values as a copy of the
503  global values.
504 
505  Inside the Vector @a values, the index `i+vdim*j` corresponds to the
506  `i`-th vector component at the `j`-th quadrature point.
507  */
508  inline void GetElementValues(int idx, Vector &values) const;
509 
510  /// Return all values associated with mesh element @a idx in a DenseMatrix.
511  /** The result is stored in the DenseMatrix @a values as a reference to the
512  global values.
513 
514  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
515  `i`-th vector component at the `j`-th quadrature point.
516  */
517  inline void GetElementValues(int idx, DenseMatrix &values);
518 
519  /// Return all values associated with mesh element @a idx in a const DenseMatrix.
520  /** The result is stored in the DenseMatrix @a values as a copy of the
521  global values.
522 
523  Inside the DenseMatrix @a values, the `(i,j)` entry corresponds to the
524  `i`-th vector component at the `j`-th quadrature point.
525  */
526  inline void GetElementValues(int idx, DenseMatrix &values) const;
527 
528  /// Write the QuadratureFunction to the stream @a out.
529  void Save(std::ostream &out) const;
530 };
531 
532 /// Overload operator<< for std::ostream and QuadratureFunction.
533 std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
534 
535 
536 double ZZErrorEstimator(BilinearFormIntegrator &blfi,
537  GridFunction &u,
538  GridFunction &flux,
539  Vector &error_estimates,
540  Array<int> *aniso_flags = NULL,
541  int with_subdomains = 1);
542 
543 /// Compute the Lp distance between two grid functions on the given element.
544 double ComputeElementLpDistance(double p, int i,
545  GridFunction& gf1, GridFunction& gf2);
546 
547 
548 /// Class used for extruding scalar GridFunctions
550 {
551 private:
552  int n;
553  Mesh *mesh_in;
554  Coefficient &sol_in;
555 public:
557  : n(_n), mesh_in(m), sol_in(s) { }
558  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
559  virtual ~ExtrudeCoefficient() { }
560 };
561 
562 /// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
563 GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
564  GridFunction *sol, const int ny);
565 
566 
567 // Inline methods
568 
569 inline void QuadratureFunction::SetSpace(QuadratureSpace *qspace_, int vdim_)
570 {
571  if (qspace_ != qspace)
572  {
573  if (own_qspace) { delete qspace; }
574  qspace = qspace_;
575  own_qspace = false;
576  }
577  vdim = (vdim_ < 0) ? vdim : vdim_;
579 }
580 
582  double *qf_data, int vdim_)
583 {
584  if (qspace_ != qspace)
585  {
586  if (own_qspace) { delete qspace; }
587  qspace = qspace_;
588  own_qspace = false;
589  }
590  vdim = (vdim_ < 0) ? vdim : vdim_;
591  NewDataAndSize(qf_data, vdim*qspace->GetSize());
592 }
593 
594 inline void QuadratureFunction::GetElementValues(int idx, Vector &values)
595 {
596  const int s_offset = qspace->element_offsets[idx];
597  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
598  values.NewDataAndSize(data + vdim*s_offset, vdim*sl_size);
599 }
600 
601 inline void QuadratureFunction::GetElementValues(int idx, Vector &values) const
602 {
603  const int s_offset = qspace->element_offsets[idx];
604  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
605  values.SetSize(vdim*sl_size);
606  double *q = data + vdim*s_offset;
607  for (int i = 0; i<values.Size(); i++)
608  {
609  values(i) = *(q++);
610  }
611 }
612 
614 {
615  const int s_offset = qspace->element_offsets[idx];
616  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
617  values.Reset(data + vdim*s_offset, vdim, sl_size);
618 }
619 
621  DenseMatrix &values) const
622 {
623  const int s_offset = qspace->element_offsets[idx];
624  const int sl_size = qspace->element_offsets[idx+1] - s_offset;
625  values.SetSize(vdim, sl_size);
626  double *q = data + vdim*s_offset;
627  for (int j = 0; j<sl_size; j++)
628  for (int i = 0; i<vdim; i++)
629  {
630  values(i,j) = *(q++);
631  }
632 }
633 
634 } // namespace mfem
635 
636 #endif
void SaveSTLTri(std::ostream &out, double p1[], double p2[], double p3[])
Definition: gridfunc.cpp:2451
int GetVDim() const
Get the vector dimension.
Definition: gridfunc.hpp:463
virtual ~GridFunction()
Destroys grid function.
Definition: gridfunc.hpp:390
std::ostream & operator<<(std::ostream &out, const Mesh &mesh)
Definition: mesh.cpp:8569
QuadratureFunction(QuadratureSpace *qspace_, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace.
Definition: gridfunc.hpp:417
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
void GetElementAverages(GridFunction &avgs) const
Definition: gridfunc.cpp:1054
void NewDataAndSize(double *d, int s)
Set the Vector data and size, deleting the old data, if owned.
Definition: vector.hpp:108
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: gridfunc.hpp:489
Class used for extruding scalar GridFunctions.
Definition: gridfunc.hpp:549
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:122
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:263
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, double &integral)
Definition: gridfunc.cpp:1314
const IntegrationRule & GetElementIntRule(int idx) const
Get the IntegrationRule associated with mesh element idx.
Definition: fespace.hpp:563
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains)
Definition: gridfunc.cpp:2620
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition: gridfunc.cpp:400
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
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:2774
void MakeOwner(FiniteElementCollection *_fec)
Make the GridFunction the owner of &#39;fec&#39; and &#39;fes&#39;.
Definition: gridfunc.hpp:96
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1)
Definition: gridfunc.cpp:276
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
Definition: gridfunc.cpp:1508
bool own_qspace
QuadratureSpace ownership flag.
Definition: gridfunc.hpp:406
int GetSize() const
Return the total number of quadrature points.
Definition: fespace.hpp:560
int VectorDim() const
Definition: gridfunc.cpp:291
virtual ~QuadratureFunction()
Definition: gridfunc.hpp:433
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Delta function coefficient.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
void ImposeBounds(int i, const Vector &weights, const Vector &_lo, const Vector &_hi)
Definition: gridfunc.cpp:1117
int vdim
Vector dimension.
Definition: gridfunc.hpp:405
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Definition: gridfunc.cpp:888
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:109
virtual double ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:287
void SetSpace(QuadratureSpace *qspace_, int vdim_=-1)
Change the QuadratureSpace and optionally the vector dimension.
Definition: gridfunc.hpp:569
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
Definition: gridfunc.cpp:1291
bool OwnsSpace()
Get the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:470
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:344
QuadratureFunction(QuadratureSpace *qspace_, double *qf_data, int vdim_=1)
Create a QuadratureFunction based on the given QuadratureSpace, using the external data...
Definition: gridfunc.hpp:425
void GetVectorFieldNodalValues(Vector &val, int comp) const
Definition: gridfunc.cpp:757
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2349
const FiniteElementSpace * FESpace() const
Definition: gridfunc.hpp:338
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Definition: gridfunc.cpp:2764
GridFunction(FiniteElementSpace *f, double *data)
Construct a GridFunction using previously allocated array data.
Definition: gridfunc.hpp:84
void GetDerivative(int comp, int der_comp, GridFunction &der)
Definition: gridfunc.cpp:828
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1646
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:116
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain)
Definition: gridfunc.cpp:229
virtual double ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:293
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
Definition: gridfunc.cpp:1015
void Reset(double *d, int h, int w)
Change the data array and the size of the DenseMatrix.
Definition: densemat.hpp:71
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
Definition: gridfunc.cpp:688
void SetOwnsSpace(bool own)
Set the QuadratureSpace ownership flag.
Definition: gridfunc.hpp:473
virtual double ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const
Definition: gridfunc.cpp:1864
GridFunction(const GridFunction &orig)
Copy constructor.
Definition: gridfunc.hpp:71
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:199
QuadratureSpace * GetSpace() const
Get the associated QuadratureSpace.
Definition: gridfunc.hpp:436
virtual double ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:301
QuadratureFunction()
Create an empty QuadratureFunction.
Definition: gridfunc.hpp:411
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:1201
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:182
void SetVDim(int vdim_)
Set the vector dimension, updating the size by calling Vector::SetSize().
Definition: gridfunc.hpp:466
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction. If all dofs are true, then tv will be set to point to th...
Definition: gridfunc.cpp:312
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:434
FiniteElementCollection * OwnFEC()
Definition: gridfunc.hpp:98
virtual ~ExtrudeCoefficient()
Definition: gridfunc.hpp:559
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:337
virtual double ComputeLpError(const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2175
QuadratureFunction & operator=(double value)
Redefine &#39;=&#39; for QuadratureFunction = constant.
Definition: gridfunc.cpp:2586
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:105
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:22
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
FiniteElementSpace * fes
FE space on which grid function lives.
Definition: gridfunc.hpp:31
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
Definition: gridfunc.cpp:1084
QuadratureSpace * qspace
Associated QuadratureSpace.
Definition: gridfunc.hpp:404
virtual double ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:278
FiniteElementCollection * fec
Used when the grid function is read from a file.
Definition: gridfunc.hpp:34
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1040
double GetDivergence(ElementTransformation &tr) const
Definition: gridfunc.cpp:907
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:147
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:2707
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition: gridfunc.cpp:571
Class for integration point with weight.
Definition: intrules.hpp:25
void SaveSTL(std::ostream &out, int TimesToRefine=1)
Definition: gridfunc.cpp:2471
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition: gridfunc.cpp:467
void GetElementValues(int idx, Vector &values)
Return all values associated with mesh element idx in a Vector.
Definition: gridfunc.hpp:594
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1717
void GetVectorValues(ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const
Definition: gridfunc.cpp:517
double infinity()
Define a shortcut for std::numeric_limits&lt;double&gt;::infinity()
Definition: vector.hpp:42
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1377
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition: gridfunc.cpp:174
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Definition: gridfunc.cpp:730
virtual double ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.cpp:2069
ExtrudeCoefficient(Mesh *m, Coefficient &s, int _n)
Definition: gridfunc.hpp:556
GridFunction & operator=(double value)
Redefine &#39;=&#39; for GridFunction = constant.
Definition: gridfunc.cpp:2331
void GetCurl(ElementTransformation &tr, Vector &curl) const
Definition: gridfunc.cpp:943
Vector data type.
Definition: vector.hpp:48
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Definition: gridfunc.hpp:75
Class representing the storage layout of a QuadratureFunction.
Definition: fespace.hpp:531
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:86
double * data
Definition: vector.hpp:53
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
long GetSequence() const
Return update counter (see Mesh::sequence)
Definition: fespace.hpp:517
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: gridfunc.cpp:327
Class representing a function through its values (scalar or vector) at quadrature points...
Definition: gridfunc.hpp:401
void GetValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:614
void GetBdrValuesFrom(const GridFunction &orig_func)
Definition: gridfunc.cpp:651
void Save(std::ostream &out) const
Write the QuadratureFunction to the stream out.
Definition: gridfunc.cpp:2604
void SaveVTK(std::ostream &out, const std::string &field_name, int ref)
Definition: gridfunc.cpp:2374
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Definition: gridfunc.hpp:243
virtual double GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition: gridfunc.cpp:385
void GetGradient(ElementTransformation &tr, Vector &grad) const
Definition: gridfunc.cpp:996
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
Definition: gridfunc.cpp:788