SERiF 0.0.1a
3+1D Stellar Structure and Evolution
Loading...
Searching...
No Matches
PyOperator.h
Go to the documentation of this file.
1
8#pragma once
9
10#include "mfem.hpp"
11#include "pybind11/pybind11.h"
12#include "pybind11/stl.h" // Needed for vectors, maps, sets, strings
13
14namespace py = pybind11;
15
20namespace serif::pybind {
21
83 class PyOperator : public mfem::Operator {
84 public:
85 // Inherit the constructors from the base mfem::Operator class.
86 // This allows Python classes to call e.g., super().__init__(size).
87 using mfem::Operator::Operator;
88
89 // --- Trampoline declarations for all overridable virtual functions ---
90
102 void Mult(const mfem::Vector &x, mfem::Vector &y) const override;
103
115 void MultTranspose(const mfem::Vector &x, mfem::Vector &y) const override;
116
127 void AddMult(const mfem::Vector &x, mfem::Vector &y, const mfem::real_t a = 1.0) const override;
128
139 void AddMultTranspose(const mfem::Vector &x, mfem::Vector &y, const mfem::real_t a = 1.0) const override;
140
153 Operator& GetGradient(const mfem::Vector &x) const override;
154
165 void AssembleDiagonal(mfem::Vector &diag) const override;
166
178 const mfem::Operator* GetProlongation() const override;
179
192 const mfem::Operator* GetRestriction() const override;
193
207 void RecoverFEMSolution(const mfem::Vector &X, const mfem::Vector &b, mfem::Vector &x) override;
208 };
209
210} // namespace serif::pybind
Trampoline class for mfem::Operator.
Definition PyOperator.h:83
const mfem::Operator * GetRestriction() const override
Get the restriction operator.
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.
void Mult(const mfem::Vector &x, mfem::Vector &y) const override
Perform the operator action: y = A*x.
Definition PyOperator.cpp:8
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 MultTranspose(const mfem::Vector &x, mfem::Vector &y) const override
Perform the transpose operator action: y = A^T*x.
void RecoverFEMSolution(const mfem::Vector &X, const mfem::Vector &b, mfem::Vector &x) override
Recover the FEM solution.
const mfem::Operator * GetProlongation() const override
Get the prolongation operator.
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).
Contains pybind11 helper classes and trampoline classes for interfacing C++ with Python.