SERiF 0.0.1a
3+1D Stellar Structure and Evolution
|
Trampoline class for mfem::Operator. More...
#include <PyOperator.h>
Public Member Functions | |
void | Mult (const mfem::Vector &x, mfem::Vector &y) const override |
Perform the operator action: y = A*x. | |
void | MultTranspose (const mfem::Vector &x, mfem::Vector &y) const override |
Perform the transpose operator action: y = A^T*x. | |
void | AddMult (const mfem::Vector &x, mfem::Vector &y, const mfem::real_t a=1.0) const override |
Perform the action y += a*(A*x). | |
void | AddMultTranspose (const mfem::Vector &x, mfem::Vector &y, const mfem::real_t a=1.0) const override |
Perform the action y += a*(A^T*x). | |
Operator & | GetGradient (const mfem::Vector &x) const override |
Get the gradient operator (Jacobian) at a given point x. | |
void | AssembleDiagonal (mfem::Vector &diag) const override |
Assemble the diagonal of the operator. | |
const mfem::Operator * | GetProlongation () const override |
Get the prolongation operator. | |
const mfem::Operator * | GetRestriction () const override |
Get the restriction operator. | |
void | RecoverFEMSolution (const mfem::Vector &X, const mfem::Vector &b, mfem::Vector &x) override |
Recover the FEM solution. | |
Trampoline class for mfem::Operator.
This class allows Python classes to inherit from mfem::Operator and correctly override its virtual functions. When a virtual function is called from C++, the trampoline ensures the call is forwarded to the Python implementation if one exists. This is crucial for integrating Python-defined linear or non-linear operators into MFEM's C++-based solvers and algorithms.
Definition at line 83 of file PyOperator.h.
|
override |
Perform the action y += a*(A*x).
Optional override.
x | The input vector. |
y | The vector to which a*(A*x) is added. |
a | Scalar multiplier (defaults to 1.0). |
Definition at line 33 of file PyOperator.cpp.
|
override |
Perform the action y += a*(A^T*x).
Optional override.
x | The input vector. |
y | The vector to which a*(A^T*x) is added. |
a | Scalar multiplier (defaults to 1.0). |
Definition at line 44 of file PyOperator.cpp.
|
override |
Assemble the diagonal of the operator.
For discrete operators (e.g., matrices), this method should compute and store the diagonal entries of the operator in the vector diag
. Optional override.
diag | Output vector to store the diagonal entries. |
Definition at line 64 of file PyOperator.cpp.
|
override |
Get the gradient operator (Jacobian) at a given point x.
For non-linear operators, this method should return the linearization (Jacobian) of the operator at the point x
. The returned Operator is owned by this Operator and should not be deleted by the caller. Optional override.
x | The point at which to evaluate the gradient. |
Definition at line 55 of file PyOperator.cpp.
|
override |
Get the prolongation operator.
Used in multilevel methods (e.g., AMG). Returns a pointer to the prolongation operator (interpolation from a coarser level to this level). The returned Operator is typically owned by this Operator. Optional override.
Definition at line 73 of file PyOperator.cpp.
|
override |
Get the restriction operator.
Used in multilevel methods (e.g., AMG). Returns a pointer to the restriction operator (projection from this level to a coarser level). Typically, this is the transpose of the prolongation operator. The returned Operator is typically owned by this Operator. Optional override.
Definition at line 81 of file PyOperator.cpp.
|
override |
Perform the operator action: y = A*x.
This is a pure virtual function in mfem::Operator and must be overridden by any Python class inheriting from PyOperator.
x | The input vector. |
y | The output vector (result of A*x). |
Definition at line 8 of file PyOperator.cpp.
|
override |
Perform the transpose operator action: y = A^T*x.
Optional override. If not overridden, MFEM's base implementation (which may raise an error or be a no-op) will be used.
x | The input vector. |
y | The output vector (result of A^T*x). |
Definition at line 23 of file PyOperator.cpp.
|
override |
Recover the FEM solution.
For operators that are part of a system solve (e.g., static condensation), this method can be used to reconstruct the full finite element solution x
from a reduced solution X
and the right-hand side b
. Optional override.
X | The reduced solution vector. |
b | The right-hand side vector. |
x | Output vector for the full FEM solution. |
Definition at line 89 of file PyOperator.cpp.