MFEM
v3.3
Finite element discretization library
|
Class for parallel grid function. More...
#include <pgridfunc.hpp>
Public Member Functions | |
ParGridFunction () | |
ParGridFunction (ParFiniteElementSpace *pf) | |
ParGridFunction (ParFiniteElementSpace *pf, double *data) | |
Construct a ParGridFunction using previously allocated array data. More... | |
ParGridFunction (ParFiniteElementSpace *pf, GridFunction *gf) | |
Construct a ParGridFunction using a GridFunction as external data. More... | |
ParGridFunction (ParFiniteElementSpace *pf, HypreParVector *tv) | |
ParGridFunction (ParMesh *pmesh, GridFunction *gf, int *partitioning=NULL) | |
ParGridFunction & | operator= (double value) |
Assign constant values to the ParGridFunction data. More... | |
ParGridFunction & | operator= (const Vector &v) |
Copy the data from a Vector to the ParGridFunction data. More... | |
ParFiniteElementSpace * | ParFESpace () const |
void | Update () |
void | SetSpace (ParFiniteElementSpace *f) |
void | MakeRef (ParFiniteElementSpace *f, double *v) |
Make the ParGridFunction reference external data on a new ParFiniteElementSpace. More... | |
void | MakeRef (ParFiniteElementSpace *f, Vector &v, int v_offset) |
Make the ParGridFunction reference external data on a new ParFiniteElementSpace. More... | |
void | Distribute (const Vector *tv) |
void | Distribute (const Vector &tv) |
void | AddDistribute (double a, const Vector *tv) |
void | AddDistribute (double a, const Vector &tv) |
virtual void | SetFromTrueDofs (const Vector &tv) |
Set the GridFunction from the given true-dof vector. More... | |
ParGridFunction & | operator= (const HypreParVector &tv) |
Short semantic for Distribute. More... | |
HypreParVector * | GetTrueDofs () const |
Returns the true dofs in a new HypreParVector. More... | |
void | ParallelAverage (Vector &tv) const |
Returns the vector averaged on the true dofs. More... | |
void | ParallelAverage (HypreParVector &tv) const |
Returns the vector averaged on the true dofs. More... | |
HypreParVector * | ParallelAverage () const |
Returns a new vector averaged on the true dofs. More... | |
void | ParallelProject (Vector &tv) const |
Returns the vector restricted to the true dofs. More... | |
void | ParallelProject (HypreParVector &tv) const |
Returns the vector restricted to the true dofs. More... | |
HypreParVector * | ParallelProject () const |
Returns a new vector restricted to the true dofs. More... | |
void | ParallelAssemble (Vector &tv) const |
Returns the vector assembled on the true dofs. More... | |
void | ParallelAssemble (HypreParVector &tv) const |
Returns the vector assembled on the true dofs. More... | |
HypreParVector * | ParallelAssemble () const |
Returns a new vector assembled on the true dofs. More... | |
void | ExchangeFaceNbrData () |
Vector & | FaceNbrData () |
const Vector & | FaceNbrData () const |
virtual double | GetValue (int i, const IntegrationPoint &ip, int vdim=1) const |
double | GetValue (ElementTransformation &T) |
void | ProjectCoefficient (Coefficient &coeff) |
void | ProjectDiscCoefficient (VectorCoefficient &coeff) |
double | ComputeL1Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const |
double | ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeL2Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const |
double | ComputeL2Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeL2Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const |
double | ComputeMaxError (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const |
double | ComputeMaxError (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeMaxError (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeLpError (const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const |
double | ComputeLpError (const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const |
virtual void | ComputeFlux (BilinearFormIntegrator &blfi, GridFunction &flux, int wcoef=1, int subdomain=-1) |
virtual void | Save (std::ostream &out) const |
void | SaveAsOne (std::ostream &out=std::cout) |
Merge the local grid functions. More... | |
virtual | ~ParGridFunction () |
Public Member Functions inherited from mfem::GridFunction | |
GridFunction () | |
GridFunction (const GridFunction &orig) | |
Copy constructor. More... | |
GridFunction (FiniteElementSpace *f) | |
Construct a GridFunction associated with the FiniteElementSpace *f. More... | |
GridFunction (FiniteElementSpace *f, double *data) | |
Construct a GridFunction using previously allocated array data. More... | |
GridFunction (Mesh *m, std::istream &input) | |
Construct a GridFunction on the given Mesh, using the data from input. More... | |
GridFunction (Mesh *m, GridFunction *gf_array[], int num_pieces) | |
void | MakeOwner (FiniteElementCollection *_fec) |
Make the GridFunction the owner of 'fec' and 'fes'. More... | |
FiniteElementCollection * | OwnFEC () |
int | VectorDim () const |
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 the data of *this . More... | |
void | GetNodalValues (int i, Array< double > &nval, int vdim=1) const |
Returns the values in the vertices of i'th element for dimension vdim. More... | |
void | GetVectorValue (int i, const IntegrationPoint &ip, Vector &val) const |
void | GetValues (int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const |
void | GetValues (int i, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const |
int | GetFaceValues (int i, int side, const IntegrationRule &ir, Vector &vals, DenseMatrix &tr, int vdim=1) const |
void | GetVectorValues (ElementTransformation &T, const IntegrationRule &ir, DenseMatrix &vals) const |
void | GetVectorValues (int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const |
int | GetFaceVectorValues (int i, int side, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr) const |
void | GetValuesFrom (GridFunction &) |
void | GetBdrValuesFrom (GridFunction &) |
void | GetVectorFieldValues (int i, const IntegrationRule &ir, DenseMatrix &vals, DenseMatrix &tr, int comp=0) const |
void | ReorderByNodes () |
For a vector grid function, makes sure that the ordering is byNODES. More... | |
void | GetNodalValues (Vector &nval, int vdim=1) const |
Return the values as a vector on mesh vertices for dimension vdim. More... | |
void | GetVectorFieldNodalValues (Vector &val, int comp) const |
void | ProjectVectorFieldOn (GridFunction &vec_field, int comp=0) |
void | GetDerivative (int comp, int der_comp, GridFunction &der) |
double | GetDivergence (ElementTransformation &tr) |
void | GetCurl (ElementTransformation &tr, Vector &curl) |
void | GetGradient (ElementTransformation &tr, Vector &grad) |
void | GetGradients (const int elem, const IntegrationRule &ir, DenseMatrix &grad) |
void | GetVectorGradient (ElementTransformation &tr, DenseMatrix &grad) |
void | GetElementAverages (GridFunction &avgs) |
void | ImposeBounds (int i, const Vector &weights, const Vector &_lo, const Vector &_hi) |
void | ImposeBounds (int i, const Vector &weights, double _min=0.0, double _max=std::numeric_limits< double >::infinity()) |
void | ProjectGridFunction (const GridFunction &src) |
void | ProjectCoefficient (Coefficient &coeff) |
void | ProjectCoefficient (Coefficient &coeff, Array< int > &dofs, int vd=0) |
void | ProjectCoefficient (VectorCoefficient &vcoeff) |
void | ProjectCoefficient (VectorCoefficient &vcoeff, Array< int > &dofs) |
void | ProjectCoefficient (Coefficient *coeff[]) |
void | ProjectDiscCoefficient (VectorCoefficient &coeff) |
void | ProjectBdrCoefficient (Coefficient &coeff, Array< int > &attr) |
void | ProjectBdrCoefficient (Coefficient *coeff[], Array< int > &attr) |
void | ProjectBdrCoefficientNormal (VectorCoefficient &vcoeff, Array< int > &bdr_attr) |
void | ProjectBdrCoefficientTangent (VectorCoefficient &vcoeff, Array< int > &bdr_attr) |
double | ComputeL2Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeL2Error (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const |
double | ComputeL2Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL, Array< int > *elems=NULL) const |
double | ComputeH1Error (Coefficient *exsol, VectorCoefficient *exgrad, Coefficient *ell_coef, double Nu, int norm_type) const |
double | ComputeMaxError (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeMaxError (Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const |
double | ComputeMaxError (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeL1Error (Coefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeW11Error (Coefficient *exsol, VectorCoefficient *exgrad, int norm_type, Array< int > *elems=NULL, const IntegrationRule *irs[]=NULL) const |
double | ComputeL1Error (VectorCoefficient &exsol, const IntegrationRule *irs[]=NULL) const |
double | ComputeLpError (const double p, Coefficient &exsol, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const |
double | ComputeLpError (const double p, VectorCoefficient &exsol, Coefficient *weight=NULL, VectorCoefficient *v_weight=NULL, const IntegrationRule *irs[]=NULL) const |
GridFunction & | operator= (double value) |
Redefine '=' for GridFunction = constant. More... | |
GridFunction & | operator= (const Vector &v) |
Copy the data from v. More... | |
GridFunction & | operator= (const GridFunction &v) |
Copy the data from v. More... | |
void | Update () |
Transform by the Space UpdateMatrix (e.g., on Mesh change). More... | |
FiniteElementSpace * | FESpace () |
const FiniteElementSpace * | FESpace () const |
void | SetSpace (FiniteElementSpace *f) |
void | MakeRef (FiniteElementSpace *f, double *v) |
Make the GridFunction reference external data on a new FiniteElementSpace. More... | |
void | MakeRef (FiniteElementSpace *f, Vector &v, int v_offset) |
Make the GridFunction reference external data on a new FiniteElementSpace. More... | |
void | SaveVTK (std::ostream &out, const std::string &field_name, int ref) |
void | SaveSTL (std::ostream &out, int TimesToRefine=1) |
virtual | ~GridFunction () |
Destroys grid function. More... | |
Public Member Functions inherited from mfem::Vector | |
Vector () | |
Default constructor for Vector. Sets size = 0 and data = NULL. More... | |
Vector (const Vector &) | |
Copy constructor. Allocates a new data array and copies the data. More... | |
Vector (int s) | |
Creates vector of size s. More... | |
Vector (double *_data, int _size) | |
Creates a vector referencing an array of doubles, owned by someone else. More... | |
void | Load (std::istream **in, int np, int *dim) |
Reads a vector from multiple files. More... | |
void | Load (std::istream &in, int Size) |
Load a vector from an input stream. More... | |
void | Load (std::istream &in) |
Load a vector from an input stream, reading the size from the stream. More... | |
void | SetSize (int s) |
Resize the vector to size s. More... | |
void | SetData (double *d) |
void | SetDataAndSize (double *d, int s) |
Set the Vector data and size. More... | |
void | NewDataAndSize (double *d, int s) |
Set the Vector data and size, deleting the old data, if owned. More... | |
void | MakeDataOwner () |
void | Destroy () |
Destroy a vector. More... | |
int | Size () const |
Returns the size of the vector. More... | |
int | Capacity () const |
Return the size of the currently allocated data array. More... | |
double * | GetData () const |
operator double * () | |
Conversion to double * . More... | |
operator const double * () const | |
Conversion to const double * . More... | |
bool | OwnsData () const |
void | StealData (double **p) |
Changes the ownership of the data; after the call the Vector is empty. More... | |
double * | StealData () |
Changes the ownership of the data; after the call the Vector is empty. More... | |
double & | Elem (int i) |
Access Vector entries. Index i = 0 .. size-1. More... | |
const double & | Elem (int i) const |
Read only access to Vector entries. Index i = 0 .. size-1. More... | |
double & | operator() (int i) |
Access Vector entries using () for 0-based indexing. More... | |
const double & | operator() (int i) const |
Read only access to Vector entries using () for 0-based indexing. More... | |
double | operator* (const double *) const |
Dot product with a double * array. More... | |
double | operator* (const Vector &v) const |
Return the inner-product. More... | |
Vector & | operator= (const double *v) |
Vector & | operator= (const Vector &v) |
Redefine '=' for vector = vector. More... | |
Vector & | operator= (double value) |
Redefine '=' for vector = constant. More... | |
Vector & | operator*= (double c) |
Vector & | operator/= (double c) |
Vector & | operator-= (double c) |
Vector & | operator-= (const Vector &v) |
Vector & | operator+= (const Vector &v) |
Vector & | Add (const double a, const Vector &Va) |
(*this) += a * Va More... | |
Vector & | Set (const double a, const Vector &x) |
(*this) = a * x More... | |
void | SetVector (const Vector &v, int offset) |
void | Neg () |
(*this) = -(*this) More... | |
void | Swap (Vector &other) |
Swap the contents of two Vectors. More... | |
void | median (const Vector &lo, const Vector &hi) |
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi. More... | |
void | GetSubVector (const Array< int > &dofs, Vector &elemvect) const |
void | GetSubVector (const Array< int > &dofs, double *elem_data) const |
void | SetSubVector (const Array< int > &dofs, const double value) |
Set the entries listed in dofs to the given value . More... | |
void | SetSubVector (const Array< int > &dofs, const Vector &elemvect) |
void | SetSubVector (const Array< int > &dofs, double *elem_data) |
void | AddElementVector (const Array< int > &dofs, const Vector &elemvect) |
Add (element) subvector to the vector. More... | |
void | AddElementVector (const Array< int > &dofs, double *elem_data) |
void | AddElementVector (const Array< int > &dofs, const double a, const Vector &elemvect) |
void | SetSubVectorComplement (const Array< int > &dofs, const double val) |
Set all vector entries NOT in the 'dofs' array to the given 'val'. More... | |
void | Print (std::ostream &out=std::cout, int width=8) const |
Prints vector to stream out. More... | |
void | Print_HYPRE (std::ostream &out) const |
Prints vector to stream out in HYPRE_Vector format. More... | |
void | Randomize (int seed=0) |
Set random values in the vector. More... | |
double | Norml2 () const |
Returns the l2 norm of the vector. More... | |
double | Normlinf () const |
Returns the l_infinity norm of the vector. More... | |
double | Norml1 () const |
Returns the l_1 norm of the vector. More... | |
double | Normlp (double p) const |
Returns the l_p norm of the vector. More... | |
double | Max () const |
Returns the maximal element of the vector. More... | |
double | Min () const |
Returns the minimal element of the vector. More... | |
double | Sum () const |
Return the sum of the vector entries. More... | |
double | DistanceTo (const double *p) const |
Compute the Euclidean distance to another vector. More... | |
int | CheckFinite () const |
virtual | ~Vector () |
Destroys vector. More... | |
Vector (N_Vector nv) | |
Construct a wrapper Vector from SUNDIALS N_Vector. More... | |
virtual N_Vector | ToNVector () |
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. More... | |
virtual void | ToNVector (N_Vector &nv) |
Update an existing wrapper SUNDIALS N_Vector to point to this Vector. More... | |
Protected Attributes | |
ParFiniteElementSpace * | pfes |
Vector | face_nbr_data |
Protected Attributes inherited from mfem::GridFunction | |
FiniteElementSpace * | fes |
FE space on which grid function lives. More... | |
FiniteElementCollection * | fec |
Used when the grid function is read from a file. More... | |
long | sequence |
Protected Attributes inherited from mfem::Vector | |
int | size |
int | allocsize |
double * | data |
Additional Inherited Members | |
Protected Member Functions inherited from mfem::GridFunction | |
void | SaveSTLTri (std::ostream &out, double p1[], double p2[], double p3[]) |
void | GetVectorGradientHat (ElementTransformation &T, DenseMatrix &gh) |
void | ProjectDeltaCoefficient (DeltaCoefficient &delta_coeff, double &integral) |
void | SumFluxAndCount (BilinearFormIntegrator &blfi, GridFunction &flux, Array< int > &counts, int wcoef, int subdomain) |
void | ProjectDiscCoefficient (VectorCoefficient &coeff, Array< int > &dof_attr) |
void | Destroy () |
Class for parallel grid function.
Definition at line 31 of file pgridfunc.hpp.
|
inline |
Definition at line 39 of file pgridfunc.hpp.
|
inline |
Definition at line 41 of file pgridfunc.hpp.
|
inline |
Construct a ParGridFunction using previously allocated array data.
The ParGridFunction does not assume ownership of data which is assumed to be of size at least pf->GetVSize()
. Similar to the GridFunction and Vector constructors for externally allocated array, the pointer data can be NULL. The data array can be replaced later using the method SetData().
Definition at line 50 of file pgridfunc.hpp.
mfem::ParGridFunction::ParGridFunction | ( | ParFiniteElementSpace * | pf, |
GridFunction * | gf | ||
) |
Construct a ParGridFunction using a GridFunction as external data.
The parallel space *pf and the space used by *gf should match. The data from *gf is used as the local data of the ParGridFunction on each processor. The ParGridFunction does not assume ownership of the data.
Definition at line 24 of file pgridfunc.cpp.
mfem::ParGridFunction::ParGridFunction | ( | ParFiniteElementSpace * | pf, |
HypreParVector * | tv | ||
) |
Creates grid function on (all) dofs from a given vector on the true dofs, i.e. P tv.
Definition at line 30 of file pgridfunc.cpp.
mfem::ParGridFunction::ParGridFunction | ( | ParMesh * | pmesh, |
GridFunction * | gf, | ||
int * | partitioning = NULL |
||
) |
Construct a ParGridFunction from the given serial GridFunction. If partitioning == NULL (default), the data from 'gf' is NOT copied.
Definition at line 36 of file pgridfunc.cpp.
|
inlinevirtual |
Definition at line 250 of file pgridfunc.hpp.
void mfem::ParGridFunction::AddDistribute | ( | double | a, |
const Vector * | tv | ||
) |
Definition at line 97 of file pgridfunc.cpp.
|
inline |
Definition at line 102 of file pgridfunc.hpp.
|
virtual |
Reimplemented from mfem::GridFunction.
Definition at line 497 of file pgridfunc.cpp.
|
inline |
Definition at line 162 of file pgridfunc.hpp.
|
inline |
Definition at line 169 of file pgridfunc.hpp.
|
inline |
Definition at line 173 of file pgridfunc.hpp.
|
inline |
Definition at line 177 of file pgridfunc.hpp.
|
inline |
Definition at line 184 of file pgridfunc.hpp.
|
inline |
Definition at line 188 of file pgridfunc.hpp.
|
inline |
Definition at line 218 of file pgridfunc.hpp.
|
inline |
When given a vector weight, compute the pointwise (scalar) error as the dot product of the vector error with the vector weight. Otherwise, the scalar error is the l_2 norm of the vector error.
Definition at line 229 of file pgridfunc.hpp.
|
inline |
Definition at line 196 of file pgridfunc.hpp.
|
inline |
Definition at line 204 of file pgridfunc.hpp.
|
inline |
Definition at line 211 of file pgridfunc.hpp.
void mfem::ParGridFunction::Distribute | ( | const Vector * | tv | ) |
Set the grid function on (all) dofs from a given vector on the true dofs, i.e. P tv.
Definition at line 92 of file pgridfunc.cpp.
|
inline |
Definition at line 100 of file pgridfunc.hpp.
void mfem::ParGridFunction::ExchangeFaceNbrData | ( | ) |
Definition at line 162 of file pgridfunc.cpp.
|
inline |
Definition at line 144 of file pgridfunc.hpp.
|
inline |
Definition at line 145 of file pgridfunc.hpp.
HypreParVector * mfem::ParGridFunction::GetTrueDofs | ( | ) | const |
Returns the true dofs in a new HypreParVector.
Definition at line 102 of file pgridfunc.cpp.
|
virtual |
Reimplemented from mfem::GridFunction.
Definition at line 213 of file pgridfunc.cpp.
|
inline |
Definition at line 150 of file pgridfunc.hpp.
void mfem::ParGridFunction::MakeRef | ( | ParFiniteElementSpace * | f, |
double * | v | ||
) |
Make the ParGridFunction reference external data on a new ParFiniteElementSpace.
This method changes the ParFiniteElementSpace associated with the ParGridFunction and sets the pointer v as external data in the ParGridFunction.
Definition at line 78 of file pgridfunc.cpp.
void mfem::ParGridFunction::MakeRef | ( | ParFiniteElementSpace * | f, |
Vector & | v, | ||
int | v_offset | ||
) |
Make the ParGridFunction reference external data on a new ParFiniteElementSpace.
This method changes the ParFiniteElementSpace associated with the ParGridFunction and sets the data of the Vector v (plus the v_offset) as external data in the ParGridFunction.
Definition at line 85 of file pgridfunc.cpp.
|
inline |
Assign constant values to the ParGridFunction data.
Definition at line 68 of file pgridfunc.hpp.
|
inline |
Copy the data from a Vector to the ParGridFunction data.
Definition at line 72 of file pgridfunc.hpp.
|
inline |
Short semantic for Distribute.
Definition at line 108 of file pgridfunc.hpp.
void mfem::ParGridFunction::ParallelAssemble | ( | Vector & | tv | ) | const |
Returns the vector assembled on the true dofs.
Definition at line 145 of file pgridfunc.cpp.
void mfem::ParGridFunction::ParallelAssemble | ( | HypreParVector & | tv | ) | const |
Returns the vector assembled on the true dofs.
Definition at line 150 of file pgridfunc.cpp.
HypreParVector * mfem::ParGridFunction::ParallelAssemble | ( | ) | const |
Returns a new vector assembled on the true dofs.
Definition at line 155 of file pgridfunc.cpp.
void mfem::ParGridFunction::ParallelAverage | ( | Vector & | tv | ) | const |
Returns the vector averaged on the true dofs.
Definition at line 109 of file pgridfunc.cpp.
void mfem::ParGridFunction::ParallelAverage | ( | HypreParVector & | tv | ) | const |
Returns the vector averaged on the true dofs.
Definition at line 115 of file pgridfunc.cpp.
HypreParVector * mfem::ParGridFunction::ParallelAverage | ( | ) | const |
Returns a new vector averaged on the true dofs.
Definition at line 121 of file pgridfunc.cpp.
void mfem::ParGridFunction::ParallelProject | ( | Vector & | tv | ) | const |
Returns the vector restricted to the true dofs.
Definition at line 128 of file pgridfunc.cpp.
void mfem::ParGridFunction::ParallelProject | ( | HypreParVector & | tv | ) | const |
Returns the vector restricted to the true dofs.
Definition at line 133 of file pgridfunc.cpp.
HypreParVector * mfem::ParGridFunction::ParallelProject | ( | ) | const |
Returns a new vector restricted to the true dofs.
Definition at line 138 of file pgridfunc.cpp.
|
inline |
Definition at line 75 of file pgridfunc.hpp.
void mfem::ParGridFunction::ProjectCoefficient | ( | Coefficient & | coeff | ) |
Definition at line 249 of file pgridfunc.cpp.
void mfem::ParGridFunction::ProjectDiscCoefficient | ( | VectorCoefficient & | coeff | ) |
Project a discontinuous vector coefficient as a grid function on a continuous parallel finite element space. The values in shared dofs are determined from the element with maximal attribute.
Definition at line 270 of file pgridfunc.cpp.
|
virtual |
Save the local portion of the ParGridFunction. It differs from the serial GridFunction::Save in that it takes into account the signs of the local dofs.
Reimplemented from mfem::GridFunction.
Definition at line 314 of file pgridfunc.cpp.
void mfem::ParGridFunction::SaveAsOne | ( | std::ostream & | out = std::cout | ) |
Merge the local grid functions.
Definition at line 329 of file pgridfunc.cpp.
|
inlinevirtual |
Set the GridFunction from the given true-dof vector.
Reimplemented from mfem::GridFunction.
Definition at line 105 of file pgridfunc.hpp.
void mfem::ParGridFunction::SetSpace | ( | ParFiniteElementSpace * | f | ) |
Definition at line 71 of file pgridfunc.cpp.
void mfem::ParGridFunction::Update | ( | ) |
Definition at line 65 of file pgridfunc.cpp.
|
protected |
Definition at line 36 of file pgridfunc.hpp.
|
protected |
Definition at line 34 of file pgridfunc.hpp.