MFEM
v3.4
Finite element discretization library
|
Class for parallel linear form. More...
#include <plinearform.hpp>
Public Member Functions | |
ParLinearForm () | |
ParLinearForm (ParFiniteElementSpace *pf) | |
ParFiniteElementSpace * | ParFESpace () const |
void | Update (ParFiniteElementSpace *pf=NULL) |
void | Update (ParFiniteElementSpace *pf, Vector &v, int v_offset) |
void | ParallelAssemble (Vector &tv) |
Assemble the vector on the true dofs, i.e. P^t v. More... | |
HypreParVector * | ParallelAssemble () |
Returns the vector assembled on the true dofs, i.e. P^t v. More... | |
double | operator() (const ParGridFunction &gf) const |
Return the action of the ParLinearForm as a linear mapping. More... | |
Public Member Functions inherited from mfem::LinearForm | |
LinearForm (FiniteElementSpace *f) | |
Creates linear form associated with FE space *f. More... | |
LinearForm () | |
FiniteElementSpace * | GetFES () |
(DEPRECATED) Return the FE space associated with the LinearForm. More... | |
FiniteElementSpace * | FESpace () |
Read+write access to the associated FiniteElementSpace. More... | |
const FiniteElementSpace * | FESpace () const |
Read-only access to the associated FiniteElementSpace. More... | |
void | AddDomainIntegrator (LinearFormIntegrator *lfi) |
Adds new Domain Integrator. More... | |
void | AddBoundaryIntegrator (LinearFormIntegrator *lfi) |
Adds new Boundary Integrator. More... | |
void | AddBdrFaceIntegrator (LinearFormIntegrator *lfi) |
Adds new Boundary Face Integrator. More... | |
void | AddBdrFaceIntegrator (LinearFormIntegrator *lfi, Array< int > &bdr_attr_marker) |
Add new Boundary Face Integrator, restricted to the given boundary attributes. More... | |
void | Assemble () |
Assembles the linear form i.e. sums over all domain/bdr integrators. More... | |
void | AssembleDelta () |
Assembles delta functions of the linear form. More... | |
void | Update () |
void | Update (FiniteElementSpace *f) |
void | Update (FiniteElementSpace *f, Vector &v, int v_offset) |
double | operator() (const GridFunction &gf) const |
Return the action of the LinearForm as a linear mapping. More... | |
~LinearForm () | |
Destroys linear form. 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 |
Return a pointer to the beginning of the Vector data. More... | |
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) |
Copy Size() entries from v. More... | |
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=mfem::out, 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 | DistanceSquaredTo (const double *p) const |
Compute the square of the Euclidean distance to another vector. 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 |
Protected Attributes inherited from mfem::Vector | |
int | size |
int | allocsize |
double * | data |
Class for parallel linear form.
Definition at line 26 of file plinearform.hpp.
|
inline |
Definition at line 32 of file plinearform.hpp.
|
inline |
Definition at line 34 of file plinearform.hpp.
|
inline |
Return the action of the ParLinearForm as a linear mapping.
Linear forms are linear functionals which map ParGridFunctions to the real numbers. This method performs this mapping which in this case is equivalent as an inner product of the ParLinearForm and ParGridFunction.
Definition at line 53 of file plinearform.hpp.
void mfem::ParLinearForm::ParallelAssemble | ( | Vector & | tv | ) |
Assemble the vector on the true dofs, i.e. P^t v.
Definition at line 34 of file plinearform.cpp.
HypreParVector * mfem::ParLinearForm::ParallelAssemble | ( | ) |
Returns the vector assembled on the true dofs, i.e. P^t v.
Definition at line 39 of file plinearform.cpp.
|
inline |
Definition at line 36 of file plinearform.hpp.
void mfem::ParLinearForm::Update | ( | ParFiniteElementSpace * | pf = NULL | ) |
Definition at line 21 of file plinearform.cpp.
void mfem::ParLinearForm::Update | ( | ParFiniteElementSpace * | pf, |
Vector & | v, | ||
int | v_offset | ||
) |
Definition at line 28 of file plinearform.cpp.
|
protected |
Definition at line 29 of file plinearform.hpp.