MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
pgridfunc.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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_PGRIDFUNC
13#define MFEM_PGRIDFUNC
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_MPI
18
20#include "pfespace.hpp"
21#include "gridfunc.hpp"
22#include <iostream>
23#include <limits>
24
25namespace mfem
26{
27
28/// @brief Compute a global Lp norm from the local Lp norms computed by each
29/// processor
30///
31/// @param[in] p Real value indicating the exponent of the $L^p$ norm.
32/// To avoid domain errors p should have a positive
33/// value, either finite or infinite.
34/// @param[in] loc_norm Local $L^p$ norm as computed separately on each
35/// processor.
36/// @param[in] comm MPI Communicator
37///
38/// @return Global $L^p$ norm, returned on every processor
39///
40/// @note Quadratures with negative weights (as in some simplex integration
41/// rules in MFEM) can produce negative integrals even with
42/// non-negative integrands. To avoid returning negative norms this
43/// function uses the absolute values of the local norms.
44/// This may lead to results which are not entirely consistent with
45/// such integration rules.
46real_t GlobalLpNorm(const real_t p, real_t loc_norm, MPI_Comm comm);
47
48/// Class for parallel grid function
50{
51protected:
52 ParFiniteElementSpace *pfes; ///< Points to the same object as #fes
53
54 /** @brief Vector used to store data from face-neighbor processors,
55 initialized by ExchangeFaceNbrData(). */
57
58 /** @brief Vector used as an MPI buffer to send face-neighbor data
59 in ExchangeFaceNbrData() to neighboring processors. */
60 //TODO: Use temporary memory to avoid CUDA malloc allocation cost.
62
64 const Array<int> &attr);
65
66public:
67 ParGridFunction() { pfes = NULL; }
68
69 /// Copy constructor. The internal vector #face_nbr_data is not copied.
71 : GridFunction(orig), pfes(orig.pfes) { }
72
74
75 /// Construct a ParGridFunction using previously allocated array @a data.
76 /** The ParGridFunction does not assume ownership of @a data which is assumed
77 to be of size at least `pf->GetVSize()`. Similar to the GridFunction and
78 Vector constructors for externally allocated array, the pointer @a data
79 can be NULL. The data array can be replaced later using the method
80 SetData().
81 */
84
85 /** @brief Construct a ParGridFunction using previously allocated Vector
86 @a base starting at the given offset, @a base_offset. */
87 ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset = 0)
88 : GridFunction(pf, base, base_offset), pfes(pf) { }
89
90 /// Construct a ParGridFunction using a GridFunction as external data.
91 /** The parallel space @a *pf and the space used by @a *gf should match. The
92 data from @a *gf is used as the local data of the ParGridFunction on each
93 processor. The ParGridFunction does not assume ownership of the data. */
95
96 /** @brief Creates grid function on (all) dofs from a given vector on the
97 true dofs, i.e. P tv. */
99
100 /** @brief Construct a local ParGridFunction from the given *global*
101 GridFunction. If @a partitioning is NULL (default), the data from @a gf
102 is NOT copied. */
103 ParGridFunction(ParMesh *pmesh, const GridFunction *gf,
104 const int *partitioning = NULL);
105
106 /** @brief Construct a ParGridFunction on a given ParMesh, @a pmesh, reading
107 from an std::istream.
108
109 In the process, a ParFiniteElementSpace and a FiniteElementCollection are
110 constructed. The new ParGridFunction assumes ownership of both. */
111 ParGridFunction(ParMesh *pmesh, std::istream &input);
112
113 /// Copy assignment. Only the data of the base class Vector is copied.
114 /** It is assumed that this object and @a rhs use ParFiniteElementSpace%s
115 that have the same size.
116
117 @note Defining this method overwrites the implicitly defined copy
118 assignment operator. */
120 { return operator=((const Vector &)rhs); }
121
122 /// Assign constant values to the ParGridFunction data.
124 { GridFunction::operator=(value); return *this; }
125
126 /// Copy the data from a Vector to the ParGridFunction data.
128 { GridFunction::operator=(v); return *this; }
129
131
132 void Update() override;
133
134 /// Associate a new FiniteElementSpace with the ParGridFunction.
135 /** The ParGridFunction is resized using the SetSize() method. The new space
136 @a f is expected to be a ParFiniteElementSpace. */
137 void SetSpace(FiniteElementSpace *f) override;
138
139 /// Associate a new parallel space with the ParGridFunction.
141
143
144 /** @brief Make the ParGridFunction reference external data on a new
145 FiniteElementSpace. */
146 /** This method changes the FiniteElementSpace associated with the
147 ParGridFunction and sets the pointer @a v as external data in the
148 ParGridFunction. The new space @a f is expected to be a
149 ParFiniteElementSpace. */
150 void MakeRef(FiniteElementSpace *f, real_t *v) override;
151
152 /** @brief Make the ParGridFunction reference external data on a new
153 ParFiniteElementSpace. */
154 /** This method changes the ParFiniteElementSpace associated with the
155 ParGridFunction and sets the pointer @a v as external data in the
156 ParGridFunction. */
158
159 /** @brief Make the ParGridFunction reference external data on a new
160 FiniteElementSpace. */
161 /** This method changes the FiniteElementSpace associated with the
162 ParGridFunction and sets the data of the Vector @a v (plus the @a
163 v_offset) as external data in the ParGridFunction. The new space @a f is
164 expected to be a ParFiniteElementSpace.
165 @note This version of the method will also perform bounds checks when
166 the build option MFEM_DEBUG is enabled. */
167 void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset) override;
168
169 /** @brief Make the ParGridFunction reference external data on a new
170 ParFiniteElementSpace. */
171 /** This method changes the ParFiniteElementSpace associated with the
172 ParGridFunction and sets the data of the Vector @a v (plus the
173 @a v_offset) as external data in the ParGridFunction.
174 @note This version of the method will also perform bounds checks when
175 the build option MFEM_DEBUG is enabled. */
176 void MakeRef(ParFiniteElementSpace *f, Vector &v, int v_offset);
177
178 /** Set the grid function on (all) dofs from a given vector on the
179 true dofs, i.e. P tv. */
180 void Distribute(const Vector *tv);
181 void Distribute(const Vector &tv) { Distribute(&tv); }
182 void AddDistribute(real_t a, const Vector *tv);
183 void AddDistribute(real_t a, const Vector &tv) { AddDistribute(a, &tv); }
184
185 /// Set the GridFunction from the given true-dof vector.
186 void SetFromTrueDofs(const Vector &tv) override { Distribute(tv); }
187
188 /// Short semantic for Distribute()
190 { Distribute(&tv); return (*this); }
191
193
194 /// Returns the true dofs in a new HypreParVector
196
197 /// Returns the vector averaged on the true dofs.
198 void ParallelAverage(Vector &tv) const;
199
200 /// Returns the vector averaged on the true dofs.
201 void ParallelAverage(HypreParVector &tv) const;
202
203 /// Returns a new vector averaged on the true dofs.
205
206 /// Returns the vector restricted to the true dofs.
207 void ParallelProject(Vector &tv) const;
208
209 /// Returns the vector restricted to the true dofs.
210 void ParallelProject(HypreParVector &tv) const;
211
212 /// Returns a new vector restricted to the true dofs.
214
215 /// Returns the vector assembled on the true dofs.
216 void ParallelAssemble(Vector &tv) const;
217
218 /// Returns the vector assembled on the true dofs.
219 void ParallelAssemble(HypreParVector &tv) const;
220
221 /// Returns a new vector assembled on the true dofs.
223
224 void ExchangeFaceNbrData();
226 const Vector &FaceNbrData() const { return face_nbr_data; }
227
228 // Redefine to handle the case when i is a face-neighbor element
229 real_t GetValue(int i, const IntegrationPoint &ip,
230 int vdim = 1) const override;
233
234 // Redefine to handle the case when T describes a face-neighbor element
236 int comp = 0, Vector *tr = NULL) const override;
237
238 void GetVectorValue(int i, const IntegrationPoint &ip,
239 Vector &val) const override;
240
241 // Redefine to handle the case when T describes a face-neighbor element
243 const IntegrationPoint &ip,
244 Vector &val, Vector *tr = NULL) const override;
245
246 /** @brief For each vdof, counts how many elements contain the vdof,
247 as containment is determined by FiniteElementSpace::GetElementVDofs(). */
248 void CountElementsPerVDof(Array<int> &elem_per_vdof) const override;
249
250 /// Parallel version of GridFunction::GetDerivative(); see its documentation.
251 void GetDerivative(int comp, int der_comp, ParGridFunction &der) const;
252
253 /** Sets the output vector @a dof_vals to the values of the degrees of
254 freedom of element @a el. If @a el is greater than or equal to the number
255 of local elements, it will be interpreted as a shifted index of a face
256 neighbor element. */
257 void GetElementDofValues(int el, Vector &dof_vals) const override;
258
260 void ProjectCoefficient(Coefficient &coeff) override;
261
263 /** @brief Project a discontinuous vector coefficient as a grid function on
264 a continuous finite element space. The values in shared dofs are
265 determined from the element with maximal attribute. */
266 void ProjectDiscCoefficient(VectorCoefficient &coeff) override;
267
268 void ProjectDiscCoefficient(Coefficient &coeff, AvgType type) override;
269
270 void ProjectDiscCoefficient(VectorCoefficient &vcoeff, AvgType type) override;
271
273
275 const Array<int> &attr) override
276 { ProjectBdrCoefficient(NULL, &vcoeff, attr); }
277
279 const Array<int> &attr) override
280 { ProjectBdrCoefficient(coeff, NULL, attr); }
281
283 const Array<int> &bdr_attr) override;
284
285 /// @brief Returns ||u_ex - u_h||_L1 in parallel for H1 or L2 elements
286 ///
287 /// @see GridFunction::ComputeL1Error(Coefficient *exsol[],
288 /// const IntegrationRule *irs[]) const
289 /// for more detailed documentation.
290 ///
291 /// @warning While this function is nominally equivalent to ComputeLpError,
292 /// with appropriate arguments, the returned errors may differ
293 /// noticeably because ComputeLpError uses a higher order
294 /// integration rule by default.
295 ///
296 /// @deprecated See @ref ComputeL1Error(Coefficient &exsol,
297 /// const IntegrationRule *irs[]) const
298 /// for the preferred implementation.
299 MFEM_DEPRECATED
301 const IntegrationRule *irs[] = NULL) const override
302 {
303#pragma GCC diagnostic push
304#pragma GCC diagnostic ignored "-Wdeprecated-declarations"
305 real_t glb_err = GlobalLpNorm(1.0,
307 pfes->GetComm());
308#pragma GCC diagnostic pop
309 return glb_err;
310 }
311
312 /// @brief Returns ||u_ex - u_h||_L1 in parallel for scalar fields
313 ///
314 /// @see GridFunction::ComputeL1Error(Coefficient &exsol,
315 /// const IntegrationRule *irs[]) const
316 /// for more detailed documentation.
318 const IntegrationRule *irs[] = NULL) const override
319 { return ComputeLpError(1.0, exsol, NULL, irs); }
320
321 /// @brief Returns ||u_ex - u_h||_L1 in parallel for vector fields
322 ///
323 /// @see GridFunction::ComputeL1Error(VectorCoefficient &exsol,
324 /// const IntegrationRule *irs[]) const
325 /// for more detailed documentation.
327 const IntegrationRule *irs[] = NULL) const override
328 { return ComputeLpError(1.0, exsol, NULL, NULL, irs); }
329
330 /// @brief Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements
331 ///
332 /// @see GridFunction::ComputeL2Error(Coefficient *exsol[],
333 /// const IntegrationRule *irs[],
334 /// const Array<int> *elems) const
335 /// for more detailed documentation.
337 const IntegrationRule *irs[] = NULL,
338 const Array<int> *elems = NULL) const override
339 {
340 return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
341 pfes->GetComm());
342 }
343
344 /// @brief Returns ||u_ex - u_h||_L2 in parallel for scalar fields
345 ///
346 /// @see GridFunction::ComputeL2Error(Coefficient &exsol,
347 /// const IntegrationRule *irs[],
348 /// const Array<int> *elems) const
349 /// for more detailed documentation.
351 const IntegrationRule *irs[] = NULL,
352 const Array<int> *elems = NULL) const override
353 {
354 return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
355 pfes->GetComm());
356 }
357
358 /// @brief Returns ||u_ex - u_h||_L2 in parallel for vector fields
359 ///
360 /// @see GridFunction::ComputeL2Error(VectorCoefficient &exsol,
361 /// const IntegrationRule *irs[],
362 /// const Array<int> *elems) const
363 /// for more detailed documentation.
365 const IntegrationRule *irs[] = NULL,
366 const Array<int> *elems = NULL) const override
367 {
368 return GlobalLpNorm(2.0, GridFunction::ComputeL2Error(exsol, irs, elems),
369 pfes->GetComm());
370 }
371
372 /// @brief Returns ||grad u_ex - grad u_h||_L2 in parallel for H1 or L2
373 /// elements
374 ///
375 /// @see GridFunction::ComputeGradError(VectorCoefficient *exgrad,
376 /// const IntegrationRule *irs[]) const
377 /// for more detailed documentation.
379 const IntegrationRule *irs[] = NULL) const override
380 {
381 return GlobalLpNorm(2.0, GridFunction::ComputeGradError(exgrad,irs),
382 pfes->GetComm());
383 }
384
385 /// @brief Returns ||curl u_ex - curl u_h||_L2 in parallel for ND elements
386 ///
387 /// @see GridFunction::ComputeCurlError(VectorCoefficient *excurl,
388 /// const IntegrationRule *irs[]) const
389 /// for more detailed documentation.
391 const IntegrationRule *irs[] = NULL) const override
392 {
393 return GlobalLpNorm(2.0, GridFunction::ComputeCurlError(excurl,irs),
394 pfes->GetComm());
395 }
396
397 /// @brief Returns ||div u_ex - div u_h||_L2 in parallel for RT elements
398 ///
399 /// @see GridFunction::ComputeDivError(Coefficient *exdiv,
400 /// const IntegrationRule *irs[]) const
401 /// for more detailed documentation.
403 const IntegrationRule *irs[] = NULL) const override
404 {
405 return GlobalLpNorm(2.0, GridFunction::ComputeDivError(exdiv,irs),
406 pfes->GetComm());
407 }
408
409 /// @brief Returns the Face Jumps error for L2 elements.
410 ///
411 /// Computes:
412 /// $$\sqrt{\sum_{faces} \int_f js[f] ell[f] (2 u_{ex} - u_1 - u_2)^2}$$
413 ///
414 /// Where js[f] is the jump_scaling evaluated on the face f and ell is the
415 /// average of ell_coef evaluated in the two elements sharing the face f.
416 ///
417 /// @param[in] exsol Pointer to a Coefficient object reproducing the
418 /// anticipated values of the scalar field, u_ex.
419 /// @param[in] ell_coeff Pointer to a Coefficient object used to compute
420 /// the averaged value ell in the above integral.
421 /// @param[in] jump_scaling Can be configured to provide scaling by
422 /// nu, nu/h, or nu*p^2/h
423 /// @param[in] irs Optional pointer to a custom integration rule
424 /// e.g. higher order than the default rule.
425 ///
426 /// @note Quadratures with negative weights (as in some simplex integration
427 /// rules in MFEM) can produce negative integrals even with
428 /// non-negative integrands. To avoid returning negative errors this
429 /// function uses the absolute values of the element-wise integrals.
430 /// This may lead to results which are not entirely consistent with
431 /// such integration rules.
433 Coefficient *ell_coeff,
434 JumpScaling jump_scaling,
435 const IntegrationRule *irs[]=NULL
436 ) const override;
437
438 /// Returns either the H1-seminorm or the DG Face Jumps error or both
439 /// depending on norm_type = 1, 2, 3
440 ///
441 /// @see GridFunction::ComputeH1Error(Coefficient *exsol,
442 /// VectorCoefficient *exgrad,
443 /// Coefficient *ell_coef,
444 /// real_t NU,
445 /// int norm_type) const
446 /// for more detailed documentation.
448 Coefficient *ell_coef, real_t Nu,
449 int norm_type) const override
450 {
451 return GlobalLpNorm(2.0,
452 GridFunction::ComputeH1Error(exsol,exgrad,ell_coef,
453 Nu, norm_type),
454 pfes->GetComm());
455 }
456
457 /// @brief Returns the error measured in H1-norm in parallel for H1 or L2
458 /// elements
459 ///
460 /// @see GridFunction::ComputeH1Error(Coefficient *exsol,
461 /// VectorCoefficient *exgrad,
462 /// const IntegrationRule *irs[]) const
463 /// for more detailed documentation.
465 const IntegrationRule *irs[] = NULL) const override
466 {
467 return GlobalLpNorm(2.0, GridFunction::ComputeH1Error(exsol,exgrad,irs),
468 pfes->GetComm());
469 }
470
471 /// @brief Returns the error measured H(div)-norm in parallel for RT elements
472 ///
473 /// @see GridFunction::ComputeHDivError(VectorCoefficient *exsol,
474 /// Coefficient *exdiv,
475 /// const IntegrationRule *irs[]) const
476 /// for more detailed documentation.
478 Coefficient *exdiv,
479 const IntegrationRule *irs[] = NULL) const override
480 {
481 return GlobalLpNorm(2.0, GridFunction::ComputeHDivError(exsol,exdiv,irs),
482 pfes->GetComm());
483 }
484
485 /// @brief Returns the error measured H(curl)-norm in parallel for ND
486 /// elements
487 ///
488 /// @see GridFunction::ComputeHCurlError(VectorCoefficient *exsol,
489 /// VectorCoefficient *excurl,
490 /// const IntegrationRule *irs[]) const
491 /// for more detailed documentation.
493 VectorCoefficient *excurl,
494 const IntegrationRule *irs[] = NULL) const override
495 {
496 return GlobalLpNorm(2.0,
497 GridFunction::ComputeHCurlError(exsol,excurl,irs),
498 pfes->GetComm());
499 }
500
501 /// @brief Returns Max|u_ex - u_h| error in parallel for scalar or vector
502 /// fields
503 ///
504 /// @note This implementation of the max error of a vector field computes
505 /// the max norm over vector components rather than the magnitude of
506 /// the vector.
507 ///
508 /// @see GridFunction::ComputeMaxError(Coefficient *exsol[],
509 /// const IntegrationRule *irs[]) const
510 /// for more detailed documentation.
512 const IntegrationRule *irs[] = NULL) const override
513 {
514 return GlobalLpNorm(infinity(),
516 pfes->GetComm());
517 }
518
519 /// @brief Returns Max|u_ex - u_h| error in parallel for scalar fields
520 ///
521 /// @see GridFunction::ComputeMaxError(Coefficient &exsol,
522 /// const IntegrationRule *irs[]) const
523 /// for more detailed documentation.
525 const IntegrationRule *irs[] = NULL) const override
526 {
527 return ComputeLpError(infinity(), exsol, NULL, irs);
528 }
529
530 /// @brief Returns Max|u_ex - u_h| error in parallel for vector fields
531 ///
532 /// @see GridFunction::ComputeMaxError(VectorCoefficient &exsol,
533 /// const IntegrationRule *irs[]) const
534 /// for more detailed documentation.
536 const IntegrationRule *irs[] = NULL) const override
537 {
538 return ComputeLpError(infinity(), exsol, NULL, NULL, irs);
539 }
540
541 /// @brief Returns ||u_ex - u_h||_Lp in parallel for scalar fields
542 ///
543 /// @see GridFunction::ComputeLpError(const real_t p,
544 /// Coefficient &exsol,
545 /// Coefficient *weight,
546 /// const IntegrationRule *irs[],
547 /// const Array<int> *elems) const
548 /// for more detailed documentation.
550 Coefficient *weight = NULL,
551 const IntegrationRule *irs[] = NULL,
552 const Array<int> *elems = NULL) const override
553 {
554 return GlobalLpNorm(p, GridFunction::ComputeLpError(p, exsol, weight,
555 irs, elems),
556 pfes->GetComm());
557 }
558
559 /// @brief Returns ||u_ex - u_h||_Lp in parallel for vector fields
560 ///
561 /// @see GridFunction::ComputeLpError(const real_t p,
562 /// VectorCoefficient &exsol,
563 /// Coefficient *weight,
564 /// VectorCoefficient *v_weight,
565 /// const IntegrationRule *irs[]) const
566 /// for more detailed documentation.
568 Coefficient *weight = NULL,
569 VectorCoefficient *v_weight = NULL,
570 const IntegrationRule *irs[] = NULL) const override
571 {
573 p, exsol, weight, v_weight, irs), pfes->GetComm());
574 }
575
577 GridFunction &flux,
578 bool wcoef = true, int subdomain = -1) override;
579
580 /** Save the local portion of the ParGridFunction. This differs from the
581 serial GridFunction::Save in that it takes into account the signs of
582 the local dofs. */
583 void Save(std::ostream &out) const override;
584
585 /// Save the ParGridFunction to a single file (written using MPI rank 0). The
586 /// given @a precision will be used for ASCII output.
587 void SaveAsOne(const char *fname, int precision=16) const;
588
589 /// Save the ParGridFunction to files (one for each MPI rank). The files will
590 /// be given suffixes according to the MPI rank. The given @a precision will
591 /// be used for ASCII output.
592 void Save(const char *fname, int precision=16) const override;
593
594 /// @brief Returns a GridFunction on MPI rank @a save_rank that does not have
595 /// any duplication of vertices/nodes at processor boundaries.
596 ///
597 /// The @a serial_mesh is obtained using ParMesh::GetSerialMesh. Note that
598 /// the @a save_rank must be the same as that used in ParMesh::GetSerialMesh.
599 ///
600 /// @note The returned GridFunction will own the newly created
601 /// FiniteElementCollection and FiniteElementSpace objects.
602 GridFunction GetSerialGridFunction(int save_rank, Mesh &serial_mesh) const;
603
604 /// @brief Returns a GridFunction on MPI rank @a save_rank that does not have
605 /// any duplication of vertices/nodes at processor boundaries.
606 ///
607 /// The given @a serial_fes must be defined on the mesh returned by
608 /// ParMesh::GetSerialMesh (with @a save_rank ranks), for example using the
609 /// space belonging to the GridFunction obtained from @ref
610 /// ParGridFunction::GetSerialGridFunction(int,Mesh &) const.
611 ///
612 /// @note The returned GridFunction does not assume ownership of @a
613 /// serial_fes.
615 int save_rank, FiniteElementSpace &serial_fes) const;
616
617 /// Write the serial GridFunction a single file (written using MPI rank 0).
618 /// The given @a precision will be used for ASCII output.
619 void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const;
620
621#ifdef MFEM_USE_ADIOS2
622 /** Save the local portion of the ParGridFunction. This differs from the
623 serial GridFunction::Save in that it takes into account the signs of
624 the local dofs. */
625 void Save(
626 adios2stream &out, const std::string &variable_name,
628 override;
629#endif
630
631 /// Merge the local grid functions
632 void SaveAsOne(std::ostream &out = mfem::out) const;
633
634 /** @brief Return a GridFunction with the values of this, prolongated to the
635 maximum order of all elements in the mesh. */
636 std::unique_ptr<ParGridFunction> ProlongateToMaxOrder() const;
637
638 virtual ~ParGridFunction() = default;
639};
640
641
642/** Performs a global L2 projection (through a HypreBoomerAMG solve) of flux
643 from supplied discontinuous space into supplied smooth (continuous, or at
644 least conforming) space, and computes the Lp norms of the differences
645 between them on each element. This is one approach to handling conforming
646 and non-conforming elements in parallel. Returns the total error estimate. */
648 const ParGridFunction &x,
649 ParFiniteElementSpace &smooth_flux_fes,
650 ParFiniteElementSpace &flux_fes,
651 Vector &errors, int norm_p = 2, real_t solver_tol = 1e-12,
652 int solver_max_it = 200);
653
654}
655
656#endif // MFEM_USE_MPI
657
658#endif
Abstract base class BilinearFormIntegrator.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition eltrans.hpp:111
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
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
virtual real_t ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Returns Max|u_ex - u_h| error for H1 or L2 elements.
Definition gridfunc.hpp:869
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition gridfunc.cpp:233
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||exsol - u_h||_L2 for scalar or vector H1 or L2 elements.
virtual real_t ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Returns ||u_ex - u_h||_L1 for H1 or L2 elements.
Definition gridfunc.hpp:962
virtual real_t ComputeLpError(const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Returns ||u_ex - u_h||_Lp for H1 or L2 elements.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void ProjectDiscCoefficient(VectorCoefficient &coeff, Array< int > &dof_attr)
virtual real_t ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns ||curl u_ex - curl u_h||_L2 for ND elements.
virtual real_t ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(div)-norm for RT elements.
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.
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition gridfunc.cpp:363
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 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:481
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:219
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Mesh data type.
Definition mesh.hpp:64
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:334
Class for parallel grid function.
Definition pgridfunc.hpp:50
void CountElementsPerVDof(Array< int > &elem_per_vdof) const override
For each vdof, counts how many elements contain the vdof, as containment is determined by FiniteEleme...
void GetDerivative(int comp, int der_comp, ParGridFunction &der) const
Parallel version of GridFunction::GetDerivative(); see its documentation.
real_t ComputeMaxError(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns Max|u_ex - u_h| error in parallel for vector fields.
ParGridFunction(ParFiniteElementSpace *pf, Vector &base, int base_offset=0)
Construct a ParGridFunction using previously allocated Vector base starting at the given offset,...
Definition pgridfunc.hpp:87
HypreParVector * ParallelAverage() const
Returns a new vector averaged on the true dofs.
real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, real_t Nu, int norm_type) const override
real_t ComputeDGFaceJumpError(Coefficient *exsol, Coefficient *ell_coeff, JumpScaling jump_scaling, const IntegrationRule *irs[]=NULL) const override
Returns the Face Jumps error for L2 elements.
void Save(std::ostream &out) const override
real_t ComputeL1Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_L1 in parallel for scalar fields.
ParGridFunction(ParFiniteElementSpace *pf, real_t *data)
Construct a ParGridFunction using previously allocated array data.
Definition pgridfunc.hpp:82
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
void ProjectDiscCoefficient(VectorCoefficient &coeff) override
Project a discontinuous vector coefficient as a grid function on a continuous finite element space....
real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const override
real_t ComputeH1Error(Coefficient *exsol, VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const override
Returns the error measured in H1-norm in parallel for H1 or L2 elements.
void ProjectBdrCoefficient(VectorCoefficient &vcoeff, const Array< int > &attr) override
Project a VectorCoefficient on the GridFunction, modifying only DOFs on the boundary associated with ...
const Vector & FaceNbrData() const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
ParGridFunction(const ParGridFunction &orig)
Copy constructor. The internal vector face_nbr_data is not copied.
Definition pgridfunc.hpp:70
real_t ComputeLpError(const real_t p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_Lp in parallel for scalar fields.
void AddDistribute(real_t a, const Vector &tv)
real_t ComputeLpError(const real_t p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_Lp in parallel for vector fields.
void SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1) override
ParGridFunction & operator=(real_t value)
Assign constant values to the ParGridFunction data.
MFEM_DEPRECATED real_t ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_L1 in parallel for H1 or L2 elements.
ParGridFunction & operator=(const HypreParVector &tv)
Short semantic for Distribute()
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
real_t ComputeMaxError(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns Max|u_ex - u_h| error in parallel for scalar fields.
real_t ComputeHDivError(VectorCoefficient *exsol, Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
Returns the error measured H(div)-norm in parallel for RT elements.
void ProjectBdrCoefficient(Coefficient *coeff[], const Array< int > &attr) override
Project a set of Coefficients on the components of the GridFunction, modifying only DOFs on the bound...
ParFiniteElementSpace * pfes
Points to the same object as fes.
Definition pgridfunc.hpp:52
Vector send_data
Vector used as an MPI buffer to send face-neighbor data in ExchangeFaceNbrData() to neighboring proce...
Definition pgridfunc.hpp:61
HypreParVector * ParallelProject() const
Returns a new vector restricted to the true dofs.
HypreParVector * ParallelAssemble() const
Returns a new vector assembled on the true dofs.
void Distribute(const Vector &tv)
real_t ComputeL1Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const override
Returns ||u_ex - u_h||_L1 in parallel for vector fields.
ParFiniteElementSpace * ParFESpace() const
real_t ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for scalar fields.
virtual ~ParGridFunction()=default
void AddDistribute(real_t a, const Vector *tv)
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const override
real_t ComputeL2Error(VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for vector fields.
real_t ComputeGradError(VectorCoefficient *exgrad, const IntegrationRule *irs[]=NULL) const override
Returns ||grad u_ex - grad u_h||_L2 in parallel for H1 or L2 elements.
real_t ComputeDivError(Coefficient *exdiv, const IntegrationRule *irs[]=NULL) const override
Returns ||div u_ex - div u_h||_L2 in parallel for RT elements.
void SaveAsSerial(const char *fname, int precision=16, int save_rank=0) const
ParGridFunction(ParFiniteElementSpace *pf)
Definition pgridfunc.hpp:73
Vector face_nbr_data
Vector used to store data from face-neighbor processors, initialized by ExchangeFaceNbrData().
Definition pgridfunc.hpp:56
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void SaveAsOne(const char *fname, int precision=16) const
GridFunction GetSerialGridFunction(int save_rank, Mesh &serial_mesh) const
Returns a GridFunction on MPI rank save_rank that does not have any duplication of vertices/nodes at ...
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition pgridfunc.cpp:91
real_t ComputeMaxError(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const override
Returns Max|u_ex - u_h| error in parallel for scalar or vector fields.
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:97
real_t ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const override
Returns the error measured H(curl)-norm in parallel for ND elements.
void Distribute(const Vector *tv)
std::unique_ptr< ParGridFunction > ProlongateToMaxOrder() const
Return a GridFunction with the values of this, prolongated to the maximum order of all elements in th...
void GetElementDofValues(int el, Vector &dof_vals) const override
ParGridFunction & operator=(const Vector &v)
Copy the data from a Vector to the ParGridFunction data.
real_t ComputeCurlError(VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const override
Returns ||curl u_ex - curl u_h||_L2 in parallel for ND elements.
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, const Array< int > &bdr_attr) override
Project the tangential components of the given VectorCoefficient on the boundary. Only boundary attri...
real_t GetValue(ElementTransformation &T)
ParGridFunction & operator=(const ParGridFunction &rhs)
Copy assignment. Only the data of the base class Vector is copied.
Class for parallel meshes.
Definition pmesh.hpp:34
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:82
Memory< real_t > data
Definition vector.hpp:85
real_t a
Definition lissajous.cpp:41
real_t infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
Definition vector.hpp:47
real_t L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, real_t solver_tol, int solver_max_it)
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 GlobalLpNorm(const real_t p, real_t loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
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)