MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
gridfunc.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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
21#endif
22#include <limits>
23#include <ostream>
24#include <string>
25
26namespace mfem
27{
28
29/// Class for grid function - Vector with associated FE space.
30class GridFunction : public Vector
31{
32protected:
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, real_t p1[], real_t p2[], real_t p3[]);
50
51 // Project the delta coefficient without scaling and return the (local)
52 // integral of the projection.
54 real_t &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. */
67
68 /// Loading helper.
69 void LegacyNCReorder();
70
71 void Destroy();
72
73public:
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<real_t> &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 real_t 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. */
221 int comp = 0, Vector *tr = NULL) const;
222
223 /** Return a vector value from within the element indicated by the
224 ElementTransformation Object. */
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
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. */
353
354 /** @brief Compute the vector gradient with respect to the reference element
355 variable. */
357
358 /** Compute $ (\int_{\Omega} (*this) \psi_i)/(\int_{\Omega} \psi_i) $,
359 where $ \psi_i $ 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 real_t min_ = 0.0, real_t 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. */
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
431protected:
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. */
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,
450 const Array<int> &attr,
451 Array<int> &values_counter);
452
454 const Array<int> &bdr_attr,
455 Array<int> &values_counter);
456
457 // Complete the computation of averages; called e.g. after
458 // AccumulateAndCountZones().
459 void ComputeMeans(AvgType type, Array<int> &zones_per_vdof);
460
461public:
462 /** @brief For each vdof, counts how many elements contain the vdof,
463 as containment is determined by FiniteElementSpace::GetElementVDofs(). */
464 virtual void CountElementsPerVDof(Array<int> &elem_per_vdof) const;
465
466 /** @brief Project a Coefficient on the GridFunction, modifying only DOFs on
467 the boundary associated with the boundary attributes marked in the
468 @a attr array. */
470 {
471 Coefficient *coeff_p = &coeff;
472 ProjectBdrCoefficient(&coeff_p, attr);
473 }
474
475 /** @brief Project a VectorCoefficient on the GridFunction, modifying only
476 DOFs on the boundary associated with the boundary attributes marked in
477 the @a attr array. */
478 virtual void ProjectBdrCoefficient(VectorCoefficient &vcoeff,
479 const Array<int> &attr);
480
481 /** @brief Project a set of Coefficient%s on the components of the
482 GridFunction, modifying only DOFs on the boundary associated with the
483 boundary attributed marked in the @a attr array. */
484 /** If a Coefficient pointer in the array @a coeff is NULL, that component
485 will not be touched. */
486 virtual void ProjectBdrCoefficient(Coefficient *coeff[],
487 const Array<int> &attr);
488
489 /** Project the normal component of the given VectorCoefficient on
490 the boundary. Only boundary attributes that are marked in
491 'bdr_attr' are projected. Assumes RT-type VectorFE GridFunction. */
493 const Array<int> &bdr_attr);
494
495 /** @brief Project the tangential components of the given VectorCoefficient
496 on the boundary. Only boundary attributes that are marked in @a bdr_attr
497 are projected. Assumes ND-type VectorFE GridFunction. */
499 const Array<int> &bdr_attr);
500
501 virtual real_t ComputeL2Error(Coefficient *exsol[],
502 const IntegrationRule *irs[] = NULL,
503 const Array<int> *elems = NULL) const;
504
505 /// Returns ||grad u_ex - grad u_h||_L2 in element ielem for H1 or L2 elements
506 virtual real_t ComputeElementGradError(int ielem, VectorCoefficient *exgrad,
507 const IntegrationRule *irs[] = NULL) const;
508
509 /// Returns ||u_ex - u_h||_L2 for H1 or L2 elements
510 /* The @a elems input variable expects a list of markers:
511 an elem marker equal to 1 will compute the L2 error on that element
512 an elem marker equal to 0 will not compute the L2 error on that element */
514 const IntegrationRule *irs[] = NULL,
515 const Array<int> *elems = NULL) const
516 { return GridFunction::ComputeLpError(2.0, exsol, NULL, irs, elems); }
517
519 const IntegrationRule *irs[] = NULL,
520 const Array<int> *elems = NULL) const;
521
522 /// Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements
524 const IntegrationRule *irs[] = NULL) const;
525
526 /// Returns ||curl u_ex - curl u_h||_L2 for ND elements
528 const IntegrationRule *irs[] = NULL) const;
529
530 /// Returns ||div u_ex - div u_h||_L2 for RT elements
531 virtual real_t ComputeDivError(Coefficient *exdiv,
532 const IntegrationRule *irs[] = NULL) const;
533
534 /// Returns the Face Jumps error for L2 elements. The error can be weighted
535 /// by a constant nu, by nu/h, or nu*p^2/h, depending on the value of
536 /// @a jump_scaling.
538 Coefficient *ell_coeff,
539 class JumpScaling jump_scaling,
540 const IntegrationRule *irs[] = NULL)
541 const;
542
543 /// Returns the Face Jumps error for L2 elements, with 1/h scaling.
544 MFEM_DEPRECATED
546 Coefficient *ell_coeff,
547 real_t Nu,
548 const IntegrationRule *irs[] = NULL) const;
549
550 /** This method is kept for backward compatibility.
551
552 Returns either the H1-seminorm, or the DG face jumps error, or both
553 depending on norm_type = 1, 2, 3. Additional arguments for the DG face
554 jumps norm: ell_coeff: mesh-depended coefficient (weight) Nu: scalar
555 constant weight */
556 virtual real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
557 Coefficient *ell_coef, real_t Nu,
558 int norm_type) const;
559
560 /// Returns the error measured in H1-norm for H1 elements or in "broken"
561 /// H1-norm for L2 elements
562 virtual real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad,
563 const IntegrationRule *irs[] = NULL) const;
564
565 /// Returns the error measured in H(div)-norm for RT elements
567 Coefficient *exdiv,
568 const IntegrationRule *irs[] = NULL) const;
569
570 /// Returns the error measured in H(curl)-norm for ND elements
572 VectorCoefficient *excurl,
573 const IntegrationRule *irs[] = NULL) const;
574
576 const IntegrationRule *irs[] = NULL) const
577 {
578 return ComputeLpError(infinity(), exsol, NULL, irs);
579 }
580
581 virtual real_t ComputeMaxError(Coefficient *exsol[],
582 const IntegrationRule *irs[] = NULL) const;
583
585 const IntegrationRule *irs[] = NULL) const
586 {
587 return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
588 }
589
591 const IntegrationRule *irs[] = NULL) const
592 { return ComputeW11Error(*exsol, NULL, 1, NULL, irs); }
593
595 const IntegrationRule *irs[] = NULL) const
596 { return ComputeLpError(1.0, exsol, NULL, irs); }
597
598 virtual real_t ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad,
599 int norm_type, const Array<int> *elems = NULL,
600 const IntegrationRule *irs[] = NULL) const;
601
603 const IntegrationRule *irs[] = NULL) const
604 { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
605
606 /* The @a elems input variable expects a list of markers:
607 an elem marker equal to 1 will compute the L2 error on that element
608 an elem marker equal to 0 will not compute the L2 error on that element */
609 virtual real_t ComputeLpError(const real_t p, Coefficient &exsol,
610 Coefficient *weight = NULL,
611 const IntegrationRule *irs[] = NULL,
612 const Array<int> *elems = NULL) const;
613
614 /** Compute the Lp error in each element of the mesh and store the results in
615 the Vector @a error. The result should be of length number of elements,
616 for example an L2 GridFunction of order zero using map type VALUE. */
617 virtual void ComputeElementLpErrors(const real_t p, Coefficient &exsol,
618 Vector &error,
619 Coefficient *weight = NULL,
620 const IntegrationRule *irs[] = NULL
621 ) const;
622
624 Vector &error,
625 const IntegrationRule *irs[] = NULL
626 ) const
627 { ComputeElementLpErrors(1.0, exsol, error, NULL, irs); }
628
630 Vector &error,
631 const IntegrationRule *irs[] = NULL
632 ) const
633 { ComputeElementLpErrors(2.0, exsol, error, NULL, irs); }
634
636 Vector &error,
637 const IntegrationRule *irs[] = NULL
638 ) const
639 { ComputeElementLpErrors(infinity(), exsol, error, NULL, irs); }
640
641 /** When given a vector weight, compute the pointwise (scalar) error as the
642 dot product of the vector error with the vector weight. Otherwise, the
643 scalar error is the l_2 norm of the vector error. */
644 virtual real_t ComputeLpError(const real_t p, VectorCoefficient &exsol,
645 Coefficient *weight = NULL,
646 VectorCoefficient *v_weight = NULL,
647 const IntegrationRule *irs[] = NULL) const;
648
649 /** Compute the Lp error in each element of the mesh and store the results in
650 the Vector @ error. The result should be of length number of elements,
651 for example an L2 GridFunction of order zero using map type VALUE. */
652 virtual void ComputeElementLpErrors(const real_t p, VectorCoefficient &exsol,
653 Vector &error,
654 Coefficient *weight = NULL,
655 VectorCoefficient *v_weight = NULL,
656 const IntegrationRule *irs[] = NULL
657 ) const;
658
660 Vector &error,
661 const IntegrationRule *irs[] = NULL
662 ) const
663 { ComputeElementLpErrors(1.0, exsol, error, NULL, NULL, irs); }
664
666 Vector &error,
667 const IntegrationRule *irs[] = NULL
668 ) const
669 { ComputeElementLpErrors(2.0, exsol, error, NULL, NULL, irs); }
670
672 Vector &error,
673 const IntegrationRule *irs[] = NULL
674 ) const
675 { ComputeElementLpErrors(infinity(), exsol, error, NULL, NULL, irs); }
676
677 virtual void ComputeFlux(BilinearFormIntegrator &blfi,
678 GridFunction &flux,
679 bool wcoef = true, int subdomain = -1);
680
681 /// Redefine '=' for GridFunction = constant.
683
684 /// Copy the data from @a v.
685 /** The size of @a v must be equal to the size of the associated
686 FiniteElementSpace #fes. */
687 GridFunction &operator=(const Vector &v);
688
689 /// Transform by the Space UpdateMatrix (e.g., on Mesh change).
690 virtual void Update();
691
692 /** Return update counter, similar to Mesh::GetSequence(). Used to
693 check if it is up to date with the space. */
694 long GetSequence() const { return fes_sequence; }
695
697 const FiniteElementSpace *FESpace() const { return fes; }
698
699 /// Associate a new FiniteElementSpace with the GridFunction.
700 /** The GridFunction is resized using the SetSize() method. */
701 virtual void SetSpace(FiniteElementSpace *f);
702
703 using Vector::MakeRef;
704
705 /** @brief Make the GridFunction reference external data on a new
706 FiniteElementSpace. */
707 /** This method changes the FiniteElementSpace associated with the
708 GridFunction and sets the pointer @a v as external data in the
709 GridFunction. */
710 virtual void MakeRef(FiniteElementSpace *f, real_t *v);
711
712 /** @brief Make the GridFunction reference external data on a new
713 FiniteElementSpace. */
714 /** This method changes the FiniteElementSpace associated with the
715 GridFunction and sets the data of the Vector @a v (plus the @a v_offset)
716 as external data in the GridFunction.
717 @note This version of the method will also perform bounds checks when
718 the build option MFEM_DEBUG is enabled. */
719 virtual void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset);
720
721 /** @brief Associate a new FiniteElementSpace and new true-dof data with the
722 GridFunction. */
723 /** - If the prolongation matrix of @a f is trivial (i.e. its method
724 FiniteElementSpace::GetProlongationMatrix() returns NULL), then the
725 method MakeRef() is called with the same arguments.
726 - Otherwise, the method SetSpace() is called with argument @a f.
727 - The internal true-dof vector is set to reference @a tv. */
729
730 /** @brief Associate a new FiniteElementSpace and new true-dof data with the
731 GridFunction. */
732 /** - If the prolongation matrix of @a f is trivial (i.e. its method
733 FiniteElementSpace::GetProlongationMatrix() returns NULL), this method
734 calls MakeRef() with the same arguments.
735 - Otherwise, this method calls SetSpace() with argument @a f.
736 - The internal true-dof vector is set to reference the sub-vector of
737 @a tv starting at the offset @a tv_offset. */
738 void MakeTRef(FiniteElementSpace *f, Vector &tv, int tv_offset);
739
740 /// Save the GridFunction to an output stream.
741 virtual void Save(std::ostream &out) const;
742
743 /// Save the GridFunction to a file. The given @a precision will be used for
744 /// ASCII output.
745 virtual void Save(const char *fname, int precision=16) const;
746
747#ifdef MFEM_USE_ADIOS2
748 /// Save the GridFunction to a binary output stream using adios2 bp format.
749 virtual void Save(adios2stream &out, const std::string& variable_name,
752#endif
753
754 /** @brief Write the GridFunction in VTK format. Note that Mesh::PrintVTK
755 must be called first. The parameter ref > 0 must match the one used in
756 Mesh::PrintVTK. */
757 void SaveVTK(std::ostream &out, const std::string &field_name, int ref);
758
759 /** @brief Write the GridFunction in STL format. Note that the mesh dimension
760 must be 2 and that quad elements will be broken into two triangles.*/
761 void SaveSTL(std::ostream &out, int TimesToRefine = 1);
762
763 /// Destroys grid function.
764 virtual ~GridFunction() { Destroy(); }
765};
766
767
768/** Overload operator<< for std::ostream and GridFunction; valid also for the
769 derived class ParGridFunction */
770std::ostream &operator<<(std::ostream &out, const GridFunction &sol);
771
772/// Class used to specify how the jump terms in
773/// GridFunction::ComputeDGFaceJumpError are scaled.
775{
776public:
783private:
784 real_t nu;
785 JumpScalingType type;
786public:
788 : nu(nu_), type(type_) { }
789 real_t Eval(real_t h, int p) const
790 {
791 real_t val = nu;
792 if (type != CONSTANT) { val /= h; }
793 if (type == P_SQUARED_OVER_H) { val *= p*p; }
794 return val;
795 }
796};
797
798/// Overload operator<< for std::ostream and QuadratureFunction.
799std::ostream &operator<<(std::ostream &out, const QuadratureFunction &qf);
800
801
804 GridFunction &flux,
805 Vector &error_estimates,
806 Array<int> *aniso_flags = NULL,
807 int with_subdomains = 1,
808 bool with_coeff = false);
809
810/// Defines the global tensor product polynomial space used by NewZZErorrEstimator
811/**
812 * See BoundingBox(...) for a description of @a angle and @a midpoint
813 */
814void TensorProductLegendre(int dim, // input
815 int order, // input
816 const Vector &x_in, // input
817 const Vector &xmax, // input
818 const Vector &xmin, // input
819 Vector &poly, // output
820 real_t angle=0.0, // input (optional)
821 const Vector *midpoint=NULL); // input (optional)
822
823/// Defines the bounding box for the face patches used by NewZZErorrEstimator
824/**
825 * By default, BoundingBox(...) computes the parameters of a minimal bounding box
826 * for the given @a face_patch that is aligned with the physical (i.e. global)
827 * Cartesian axes. This means that the size of the bounding box will depend on the
828 * orientation of the patch. It is better to construct an orientation-independent box.
829 * This is implemented for 2D patches. The parameters @a angle and @a midpoint encode
830 * the necessary additional geometric information.
831 *
832 * @a iface : Index of the face that the patch corresponds to.
833 * This is used to compute @a angle and @a midpoint.
834 *
835 * @a angle : The angle the patch face makes with the x-axis.
836 * @a midpoint : The midpoint of the face.
837 */
838void BoundingBox(const Array<int> &face_patch, // input
839 FiniteElementSpace *ufes, // input
840 int order, // input
841 Vector &xmin, // output
842 Vector &xmax, // output
843 real_t &angle, // output
844 Vector &midpoint, // output
845 int iface=-1); // input (optional)
846
847/// A ``true'' ZZ error estimator that uses face-based patches for flux reconstruction.
848/**
849 * Only two-element face patches are ever used:
850 * - For conforming faces, the face patch consists of its two neighboring elements.
851 * - In the non-conforming setting, only the face patches associated to fine-scale
852 * element faces are used. These face patches always consist of two elements
853 * delivered by mesh::GetFaceElements(Face, *Elem1, *Elem2).
854 */
856 GridFunction &u, // input
857 Vector &error_estimates, // output
858 bool subdomain_reconstruction = true, // input (optional)
859 bool with_coeff = false, // input (optional)
860 real_t tichonov_coeff = 0.0); // input (optional)
861
862/// Compute the Lp distance between two grid functions on the given element.
864 GridFunction& gf1, GridFunction& gf2);
865
866
867/// Class used for extruding scalar GridFunctions
869{
870private:
871 int n;
872 Mesh *mesh_in;
873 Coefficient &sol_in;
874public:
876 : n(n_), mesh_in(m), sol_in(s) { }
877 virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip);
879};
880
881/// Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
882GridFunction *Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d,
883 GridFunction *sol, const int ny);
884
885} // namespace mfem
886
887#endif
Abstract base class BilinearFormIntegrator.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Delta function coefficient optionally multiplied by a weight coefficient and a scaled time dependent ...
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Class used for extruding scalar GridFunctions.
Definition gridfunc.hpp:869
ExtrudeCoefficient(Mesh *m, Coefficient &s, int n_)
Definition gridfunc.hpp:875
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
GridFunction(const GridFunction &orig)
Copy constructor. The internal true-dof vector t_vec is not copied.
Definition gridfunc.hpp:78
void ComputeMeans(AvgType type, Array< int > &zones_per_vdof)
void GetLaplacians(int i, const IntegrationRule &ir, Vector &laps, int vdim=1) const
Definition gridfunc.cpp:569
void AccumulateAndCountBdrTangentValues(VectorCoefficient &vcoeff, const Array< int > &bdr_attr, Array< int > &values_counter)
GridFunction(FiniteElementSpace *f, real_t *data)
Construct a GridFunction using previously allocated array data.
Definition gridfunc.hpp:92
virtual real_t 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:513
void ProjectVectorFieldOn(GridFunction &vec_field, int comp=0)
virtual real_t ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:590
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...
virtual real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const
Definition gridfunc.cpp:449
virtual void ComputeElementL1Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:623
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....
virtual real_t ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, class JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition gridfunc.cpp:521
virtual real_t ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
virtual real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const
void GetGradients(ElementTransformation &tr, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
FiniteElementCollection * fec
Used when the grid function is read from a file. It can also be set explicitly, see MakeOwner().
Definition gridfunc.hpp:40
virtual ~GridFunction()
Destroys grid function.
Definition gridfunc.hpp:764
void AccumulateAndCountBdrValues(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr, Array< int > &values_counter)
virtual void ComputeElementL2Errors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:629
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
virtual real_t ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:575
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition gridfunc.cpp:164
FiniteElementCollection * OwnFEC()
Definition gridfunc.hpp:124
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition gridfunc.cpp:203
void ImposeBounds(int i, const Vector &weights, const Vector &lo_, const Vector &hi_)
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:144
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
void MakeTRef(FiniteElementSpace *f, real_t *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition gridfunc.cpp:221
virtual real_t ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:594
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:150
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...
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
void GetVectorFieldValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const
virtual real_t ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:602
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition gridfunc.hpp:130
void GetElementAverages(GridFunction &avgs) const
virtual real_t 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.
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition gridfunc.hpp:122
Vector & GetTrueVector()
Read and write access to the (optional) internal true-dof Vector.
Definition gridfunc.hpp:137
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...
const FiniteElementSpace * FESpace() const
Definition gridfunc.hpp:697
virtual void ComputeElementLpErrors(const real_t p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition gridfunc.cpp:381
virtual void GetElementDofValues(int el, Vector &dof_vals) const
virtual real_t ComputeLpError(const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
GridFunction(FiniteElementSpace *f)
Construct a GridFunction associated with the FiniteElementSpace *f.
Definition gridfunc.hpp:83
void GetGradients(const int elem, const IntegrationRule &ir, DenseMatrix &grad) const
Extension of GetGradient(...) for a collection of IntegrationPoints.
Definition gridfunc.hpp:346
void SaveSTLTri(std::ostream &out, real_t p1[], real_t p2[], real_t p3[])
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
virtual void ComputeElementMaxErrors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:671
int GetFaceVectorValues(int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
virtual void ComputeElementL2Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:665
virtual void ComputeElementL1Errors(VectorCoefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:659
void GetValuesFrom(const GridFunction &orig_func)
void LegacyNCReorder()
Loading helper.
void ProjectDeltaCoefficient(DeltaCoefficient &delta_coeff, real_t &integral)
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Definition gridfunc.cpp:308
FiniteElementSpace * fes
FE space on which the grid function lives. Owned if fec is not NULL.
Definition gridfunc.hpp:34
virtual real_t ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
void GetBdrValuesFrom(const GridFunction &orig_func)
virtual real_t ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
int VectorDim() const
Definition gridfunc.cpp:323
long GetSequence() const
Definition gridfunc.hpp:694
virtual real_t ComputeW11Error(Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, const Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const
void GetDerivative(int comp, int der_comp, GridFunction &der)
Compute a certain derivative of a function's component. Derivatives of the function are computed at t...
int GetFaceValues(int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const
Definition gridfunc.cpp:664
int CurlDim() const
Definition gridfunc.cpp:346
virtual real_t ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:584
void SumFluxAndCount(BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, bool wcoef, int subdomain)
Definition gridfunc.cpp:251
GridFunction & operator=(const GridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Definition gridfunc.hpp:116
virtual real_t ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns ||div u_ex - div u_h||_L2 for RT elements.
virtual void ComputeElementMaxErrors(Coefficient &exsol, Vector &error, const IntegrationRule *irs[]=NULL) const
Definition gridfunc.hpp:635
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition gridfunc.cpp:366
virtual real_t ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const
Returns ||grad u_ex - grad u_h||_L2 for H1 or L2 elements.
void GetNodalValues(int i, Array< real_t > &nval, int vdim=1) const
Returns the values in the vertices of i'th element for dimension vdim.
Definition gridfunc.cpp:395
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
Definition gridfunc.cpp:476
real_t GetDivergence(ElementTransformation &tr) const
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 ...
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
void GetCurl(ElementTransformation &tr, Vector &curl) const
void GetGradient(ElementTransformation &tr, Vector &grad) const
Gradient of a scalar function at a quadrature point.
void GetVectorGradientHat(ElementTransformation &T, DenseMatrix &gh) const
Compute the vector gradient with respect to the reference element variable.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition gridfunc.cpp:195
void GetHessians(int i, const IntegrationRule &ir, DenseMatrix &hess, int vdim=1) const
Definition gridfunc.cpp:608
void GetVectorValues(int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const
Definition gridfunc.cpp:714
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh.
void GetVectorFieldNodalValues(Vector &val, int comp) const
void ProjectBdrCoefficient(Coefficient &coeff, const Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition gridfunc.hpp:469
void ReorderByNodes()
For a vector grid function, makes sure that the ordering is byNODES.
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
JumpScaling(real_t nu_=1.0, JumpScalingType type_=CONSTANT)
Definition gridfunc.hpp:787
real_t Eval(real_t h, int p) const
Definition gridfunc.hpp:789
Mesh data type.
Definition mesh.hpp:56
Represents values or vectors of values at quadrature points on a mesh.
Definition qfunction.hpp:24
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:80
Memory< real_t > data
Definition vector.hpp:83
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition vector.hpp:139
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:602
int dim
Definition ex24.cpp:53
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
void TensorProductLegendre(int dim, int order, const Vector &x_in, const Vector &xmax, const Vector &xmin, Vector &poly, real_t angle, const Vector *midpoint)
Defines the global tensor product polynomial space used by NewZZErorrEstimator.
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:45
GridFunction * Extrude1DGridFunction(Mesh *mesh, Mesh *mesh2d, GridFunction *sol, const int ny)
Extrude a scalar 1D GridFunction, after extruding the mesh with Extrude1D.
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
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
real_t ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains, bool with_coeff)
real_t LSZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, Vector &error_estimates, bool subdomain_reconstruction, bool with_coeff, real_t tichonov_coeff)
A `‘true’' ZZ error estimator that uses face-based patches for flux reconstruction.
void BoundingBox(const Array< int > &patch, FiniteElementSpace *ufes, int order, Vector &xmin, Vector &xmax, real_t &angle, Vector &midpoint, int iface)
Defines the bounding box for the face patches used by NewZZErorrEstimator.
real_t ComputeElementLpDistance(real_t p, int i, GridFunction &gf1, GridFunction &gf2)
Compute the Lp distance between two grid functions on the given element.
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
real_t p(const Vector &x, real_t t)
RefCoord s[3]