MFEM
v4.3.0
Finite element discretization library
|
A class representing a general parametric parallel block nonlinear operator defined on the Cartesian product of multiple ParFiniteElementSpaces. More...
#include <pparamnonlinearform.hpp>
Public Member Functions | |
virtual double | GetEnergy (const Vector &x) const |
Computes the energy of the system. More... | |
ParParametricBNLForm () | |
Construct an empty ParParametricBNLForm. Initialize with SetParSpaces(). More... | |
ParParametricBNLForm (Array< ParFiniteElementSpace * > &statef, Array< ParFiniteElementSpace * > ¶mf) | |
Construct a ParParametricBNLForm on the given set of parametric and state ParFiniteElementSpaces. More... | |
ParFiniteElementSpace * | ParFESpace (int k) |
Return the k-th parallel FE state space of the ParParametricBNLForm. More... | |
const ParFiniteElementSpace * | ParFESpace (int k) const |
Return the k-th parallel FE state space of the ParParametricBNLForm (const version). More... | |
ParFiniteElementSpace * | ParParamFESpace (int k) |
const ParFiniteElementSpace * | ParParamFESpace (int k) const |
Return the k-th parallel FE parameters space of the ParParametricBNLForm (const version). More... | |
void | SetParSpaces (Array< ParFiniteElementSpace * > &statef, Array< ParFiniteElementSpace * > ¶mf) |
Set the parallel FE spaces for the state and the parametric fields. After a call to SetParSpaces(), the essential b.c. and the gradient-type (if different from the default) must be set again. More... | |
virtual void | SetEssentialBC (const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs) |
Set the state essential BCs. Here, rhs is a true dof vector! More... | |
virtual void | SetParamEssentialBC (const Array< Array< int > * > &bdr_attr_is_ess, Array< Vector * > &rhs) |
Set the essential boundary conditions on the parametric fields. More... | |
virtual void | Mult (const Vector &x, Vector &y) const |
Calculates the residual for a state input given by block T-Vector. The result is Block T-Vector! The parametric fields should be set in advance by calling SetParamFields(). More... | |
virtual void | ParamMult (const Vector &x, Vector &y) const |
Calculates the product of the adjoint field and the derivative of the state residual with respect to the parametric fields. The adjoint and the state fields should be set in advance by calling SetAdjointFields() and SetStateFields(). The input and the result are block T-Vectors! More... | |
const BlockOperator & | GetLocalGradient (const Vector &x) const |
Return the local block gradient matrix for the given true-dof vector x. More... | |
virtual BlockOperator & | GetGradient (const Vector &x) const |
Return the block gradient matrix for the given true-dof vector x. More... | |
void | SetGradientType (Operator::Type tid) |
Set the operator type id for the blocks of the parallel gradient matrix/operator. The default type is Operator::Hypre_ParCSR. More... | |
virtual | ~ParParametricBNLForm () |
Destructor. More... | |
virtual void | SetStateFields (const Vector &xv) const |
Set the state fields. More... | |
virtual void | SetAdjointFields (const Vector &av) const |
Set the adjoint fields. More... | |
virtual void | SetParamFields (const Vector &dv) const |
Set the parameters/design fields. More... | |
Public Member Functions inherited from mfem::ParametricBNLForm | |
ParametricBNLForm () | |
Construct an empty BlockNonlinearForm. Initialize with SetSpaces(). More... | |
ParametricBNLForm (Array< FiniteElementSpace * > &statef, Array< FiniteElementSpace * > ¶mf) | |
Construct a BlockNonlinearForm on the given set of FiniteElementSpaces. More... | |
FiniteElementSpace * | FESpace (int k) |
Return the k-th FE space of the ParametricBNLForm. More... | |
FiniteElementSpace * | ParamFESpace (int k) |
Return the k-th parametric FE space of the ParametricBNLForm. More... | |
const FiniteElementSpace * | FESpace (int k) const |
Return the k-th FE space of the BlockNonlinearForm (const version). More... | |
const FiniteElementSpace * | ParamFESpace (int k) const |
Array < ParametricBNLFormIntegrator * > & | GetDNFI () |
Return the integrators. More... | |
void | SetSpaces (Array< FiniteElementSpace * > &statef, Array< FiniteElementSpace * > ¶mf) |
(Re)initialize the ParametricBNLForm. More... | |
const Array< int > & | GetBlockOffsets () const |
Return the regular dof offsets. More... | |
const Array< int > & | GetBlockTrueOffsets () const |
Return the true-dof offsets. More... | |
const Array< int > & | ParamGetBlockOffsets () const |
Return the regular dof offsets for the parameters. More... | |
const Array< int > & | ParamGetBlockTrueOffsets () const |
Return the true-dof offsets for the parameters. More... | |
void | AddDomainIntegrator (ParametricBNLFormIntegrator *nlfi) |
Adds new Domain Integrator. More... | |
void | AddInteriorFaceIntegrator (ParametricBNLFormIntegrator *nlfi) |
Adds new Interior Face Integrator. More... | |
void | AddBdrFaceIntegrator (ParametricBNLFormIntegrator *nlfi) |
Adds new Boundary Face Integrator. More... | |
void | AddBdrFaceIntegrator (ParametricBNLFormIntegrator *nlfi, Array< int > &bdr_marker) |
Adds new Boundary Face Integrator, restricted to specific boundary attributes. More... | |
virtual | ~ParametricBNLForm () |
Destructor. More... | |
Public Member Functions inherited from mfem::Operator | |
void | InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const |
Initializes memory for true vectors of linear system. More... | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). More... | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). More... | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). More... | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). More... | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). More... | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). More... | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. More... | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Action of the transpose operator: y=A^t(x) . The default behavior in class Operator is to generate an error. More... | |
virtual void | AssembleDiagonal (Vector &diag) const |
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed. More... | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More... | |
virtual const Operator * | GetOutputRestriction () const |
Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0) |
Form a constrained linear system using a matrix-free approach. More... | |
void | FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B) |
Form a column-constrained linear system using a matrix-free approach. More... | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem(). More... | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. More... | |
void | FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints). More... | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. More... | |
void | PrintMatlab (std::ostream &out, int n=0, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Type | GetType () const |
Return the type ID of the Operator class. More... | |
Protected Attributes | |
BlockVector | xs_true |
BlockVector | ys_true |
Array2D< OperatorHandle * > | phBlockGrad |
BlockOperator * | pBlockGrad |
Protected Attributes inherited from mfem::ParametricBNLForm | |
Array< FiniteElementSpace * > | fes |
FE spaces on which the form lives. More... | |
Array< FiniteElementSpace * > | paramfes |
FE spaces for the parametric fields. More... | |
int | paramheight |
int | paramwidth |
Array < ParametricBNLFormIntegrator * > | dnfi |
Set of Domain Integrators to be assembled (added). More... | |
Array < ParametricBNLFormIntegrator * > | fnfi |
Set of interior face Integrators to be assembled (added). More... | |
Array < ParametricBNLFormIntegrator * > | bfnfi |
Set of Boundary Face Integrators to be assembled (added). More... | |
Array< Array< int > * > | bfnfi_marker |
BlockVector | xs |
BlockVector | ys |
BlockVector | prmxs |
BlockVector | prmys |
BlockVector | xsv |
BlockVector | xdv |
BlockVector | adv |
Array2D< SparseMatrix * > | Grads |
Array2D< SparseMatrix * > | cGrads |
BlockOperator * | BlockGrad |
Array< int > | block_offsets |
Array< int > | block_trueOffsets |
Array< int > | paramblock_offsets |
Array< int > | paramblock_trueOffsets |
Array< Array< int > * > | ess_tdofs |
Array< Array< int > * > | paramess_tdofs |
Array< const Operator * > | P |
Array of pointers to the prolongation matrix of fes, may be NULL. More... | |
Array< const Operator * > | Pparam |
Array of pointers to the prolongation matrix of paramfes, may be NULL. More... | |
Array< const SparseMatrix * > | cP |
Array of results of dynamic-casting P to SparseMatrix pointer. More... | |
Array< const SparseMatrix * > | cPparam |
Array of results of dynamic-casting Pparam to SparseMatrix pointer. More... | |
bool | is_serial = true |
Indicator if the Operator is part of a parallel run. More... | |
bool | needs_prolongation = false |
Indicator if the Operator needs prolongation on assembly. More... | |
bool | prmneeds_prolongation = false |
Indicator if the Operator needs prolongation on assembly. More... | |
BlockVector | aux1 |
BlockVector | aux2 |
BlockVector | prmaux1 |
BlockVector | prmaux2 |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. More... | |
int | width |
Dimension of the input / number of columns in the matrix. More... | |
Additional Inherited Members | |
Public Types inherited from mfem::Operator | |
enum | DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP } |
Defines operator diagonal policy upon elimination of rows and/or columns. More... | |
enum | Type { ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ, PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE, PETSC_MATGENERIC, Complex_Operator, MFEM_ComplexSparseMat, Complex_Hypre_ParCSR } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Protected Member Functions inherited from mfem::ParametricBNLForm | |
const BlockVector & | Prolongate (const BlockVector &bx) const |
const BlockVector & | ParamProlongate (const BlockVector &bx) const |
double | GetEnergyBlocked (const BlockVector &bx, const BlockVector &dx) const |
void | MultBlocked (const BlockVector &bx, const BlockVector &dx, BlockVector &by) const |
void | MultParamBlocked (const BlockVector &bx, const BlockVector &ax, const BlockVector &dx, BlockVector &dy) const |
void | ComputeGradientBlocked (const BlockVector &bx, const BlockVector &dx) const |
Specialized version of GetGradient() for BlockVector. More... | |
Protected Member Functions inherited from mfem::Operator | |
void | FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout) |
see FormSystemOperator() More... | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() More... | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". More... | |
A class representing a general parametric parallel block nonlinear operator defined on the Cartesian product of multiple ParFiniteElementSpaces.
The ParParametricBNLForm takes as input, and returns as output, vectors on the true dofs.
Definition at line 29 of file pparamnonlinearform.hpp.
|
inline |
Construct an empty ParParametricBNLForm. Initialize with SetParSpaces().
Definition at line 41 of file pparamnonlinearform.hpp.
mfem::ParParametricBNLForm::ParParametricBNLForm | ( | Array< ParFiniteElementSpace * > & | statef, |
Array< ParFiniteElementSpace * > & | paramf | ||
) |
Construct a ParParametricBNLForm on the given set of parametric and state ParFiniteElementSpaces.
Definition at line 20 of file pparamnonlinearform.cpp.
|
virtual |
Destructor.
Definition at line 313 of file pparamnonlinearform.cpp.
|
virtual |
Computes the energy of the system.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 125 of file pparamnonlinearform.cpp.
|
virtual |
Return the block gradient matrix for the given true-dof vector x.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 245 of file pparamnonlinearform.cpp.
const BlockOperator & mfem::ParParametricBNLForm::GetLocalGradient | ( | const Vector & | x | ) | const |
Return the local block gradient matrix for the given true-dof vector x.
Return the local gradient matrix for the given true-dof vector x.
Definition at line 205 of file pparamnonlinearform.cpp.
Calculates the residual for a state input given by block T-Vector. The result is Block T-Vector! The parametric fields should be set in advance by calling SetParamFields().
Reimplemented from mfem::ParametricBNLForm.
Definition at line 144 of file pparamnonlinearform.cpp.
Calculates the product of the adjoint field and the derivative of the state residual with respect to the parametric fields. The adjoint and the state fields should be set in advance by calling SetAdjointFields() and SetStateFields(). The input and the result are block T-Vectors!
Block T-Vector to Block T-Vector.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 174 of file pparamnonlinearform.cpp.
ParFiniteElementSpace * mfem::ParParametricBNLForm::ParFESpace | ( | int | k | ) |
Return the k-th parallel FE state space of the ParParametricBNLForm.
Definition at line 67 of file pparamnonlinearform.cpp.
const ParFiniteElementSpace * mfem::ParParametricBNLForm::ParFESpace | ( | int | k | ) | const |
Return the k-th parallel FE state space of the ParParametricBNLForm (const version).
Definition at line 72 of file pparamnonlinearform.cpp.
ParFiniteElementSpace * mfem::ParParametricBNLForm::ParParamFESpace | ( | int | k | ) |
Return the k-th parallel FE parameters space of the ParParametricBNLForm.
Definition at line 78 of file pparamnonlinearform.cpp.
const ParFiniteElementSpace * mfem::ParParametricBNLForm::ParParamFESpace | ( | int | k | ) | const |
Return the k-th parallel FE parameters space of the ParParametricBNLForm (const version).
Definition at line 83 of file pparamnonlinearform.cpp.
|
virtual |
Set the adjoint fields.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 338 of file pparamnonlinearform.cpp.
|
virtual |
Set the state essential BCs. Here, rhs is a true dof vector!
Reimplemented from mfem::ParametricBNLForm.
Definition at line 89 of file pparamnonlinearform.cpp.
void mfem::ParParametricBNLForm::SetGradientType | ( | Operator::Type | tid | ) |
Set the operator type id for the blocks of the parallel gradient matrix/operator. The default type is Operator::Hypre_ParCSR.
Definition at line 234 of file pparamnonlinearform.cpp.
|
virtual |
Set the essential boundary conditions on the parametric fields.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 107 of file pparamnonlinearform.cpp.
|
virtual |
Set the parameters/design fields.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 349 of file pparamnonlinearform.cpp.
void mfem::ParParametricBNLForm::SetParSpaces | ( | Array< ParFiniteElementSpace * > & | statef, |
Array< ParFiniteElementSpace * > & | paramf | ||
) |
Set the parallel FE spaces for the state and the parametric fields. After a call to SetParSpaces(), the essential b.c. and the gradient-type (if different from the default) must be set again.
Definition at line 29 of file pparamnonlinearform.cpp.
|
virtual |
Set the state fields.
Reimplemented from mfem::ParametricBNLForm.
Definition at line 326 of file pparamnonlinearform.cpp.
|
mutableprotected |
Definition at line 34 of file pparamnonlinearform.hpp.
|
mutableprotected |
Definition at line 33 of file pparamnonlinearform.hpp.
|
mutableprotected |
Definition at line 32 of file pparamnonlinearform.hpp.
|
mutableprotected |
Definition at line 32 of file pparamnonlinearform.hpp.