|
| 1 | +/****************************************************************************** |
| 2 | +* SOFA, Simulation Open-Framework Architecture * |
| 3 | +* (c) 2021 INRIA, USTL, UJF, CNRS, MGH * |
| 4 | +* * |
| 5 | +* This program is free software; you can redistribute it and/or modify it * |
| 6 | +* under the terms of the GNU Lesser General Public License as published by * |
| 7 | +* the Free Software Foundation; either version 2.1 of the License, or (at * |
| 8 | +* your option) any later version. * |
| 9 | +* * |
| 10 | +* This program is distributed in the hope that it will be useful, but WITHOUT * |
| 11 | +* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or * |
| 12 | +* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License * |
| 13 | +* for more details. * |
| 14 | +* * |
| 15 | +* You should have received a copy of the GNU Lesser General Public License * |
| 16 | +* along with this program. If not, see <http://www.gnu.org/licenses/>. * |
| 17 | +******************************************************************************* |
| 18 | +* Contact information: [email protected] * |
| 19 | +******************************************************************************/ |
| 20 | +#include <SofaPython3/SofaLinearSystem/Binding_LinearSystem.h> |
| 21 | +#include <SofaPython3/SofaLinearSystem/Binding_LinearSystem_doc.h> |
| 22 | +#include <pybind11/pybind11.h> |
| 23 | + |
| 24 | +#include <SofaPython3/Sofa/Core/Binding_Base.h> |
| 25 | +#include <SofaPython3/PythonFactory.h> |
| 26 | + |
| 27 | +#include <sofa/component/linearsystem/TypedMatrixLinearSystem.h> |
| 28 | +#include <pybind11/eigen.h> |
| 29 | +#include <sofa/linearalgebra/CompressedRowSparseMatrix.h> |
| 30 | + |
| 31 | + |
| 32 | +namespace py { using namespace pybind11; } |
| 33 | + |
| 34 | +namespace sofapython3 { |
| 35 | + |
| 36 | +template<class Real> |
| 37 | +using EigenSparseMatrix = Eigen::SparseMatrix<Real, Eigen::RowMajor>; |
| 38 | + |
| 39 | +template<class Real> |
| 40 | +using EigenMatrixMap = Eigen::Map<Eigen::SparseMatrix<Real, Eigen::RowMajor> >; |
| 41 | + |
| 42 | +template<class Real> |
| 43 | +using Vector = Eigen::Matrix<Real,Eigen::Dynamic, 1>; |
| 44 | + |
| 45 | +template<class Real> |
| 46 | +using EigenVectorMap = Eigen::Map<Vector<Real>>; |
| 47 | + |
| 48 | +template<class TBlock> |
| 49 | +EigenSparseMatrix<typename sofa::linearalgebra::CompressedRowSparseMatrix<TBlock>::Real> |
| 50 | +toEigen(sofa::linearalgebra::CompressedRowSparseMatrix<TBlock>& matrix) |
| 51 | +{ |
| 52 | + using Real = typename sofa::linearalgebra::CompressedRowSparseMatrix<TBlock>::Real; |
| 53 | + if constexpr (std::is_same_v<TBlock, Real>) |
| 54 | + { |
| 55 | + matrix.compress(); |
| 56 | + return EigenMatrixMap<Real>(matrix.rows(), matrix.cols(), matrix.getColsValue().size(), |
| 57 | + (typename EigenMatrixMap<Real>::StorageIndex*)matrix.rowBegin.data(), |
| 58 | + (typename EigenMatrixMap<Real>::StorageIndex*)matrix.colsIndex.data(), |
| 59 | + matrix.colsValue.data()); |
| 60 | + } |
| 61 | + else |
| 62 | + { |
| 63 | + sofa::linearalgebra::CompressedRowSparseMatrix<typename TBlock::Real> filtered; |
| 64 | + filtered.copyNonZeros(matrix); |
| 65 | + filtered.compress(); |
| 66 | + return EigenMatrixMap<Real>(filtered.rows(), filtered.cols(), filtered.getColsValue().size(), |
| 67 | + (typename EigenMatrixMap<Real>::StorageIndex*)filtered.rowBegin.data(), |
| 68 | + (typename EigenMatrixMap<Real>::StorageIndex*)filtered.colsIndex.data(), |
| 69 | + filtered.colsValue.data()); |
| 70 | + } |
| 71 | +} |
| 72 | + |
| 73 | +template<class TBlock> |
| 74 | +void bindLinearSystems(py::module &m) |
| 75 | +{ |
| 76 | + using CRS = sofa::linearalgebra::CompressedRowSparseMatrix<TBlock>; |
| 77 | + using Real = typename CRS::Real; |
| 78 | + using CRSLinearSystem = sofa::component::linearsystem::TypedMatrixLinearSystem<CRS, sofa::linearalgebra::FullVector<Real> >; |
| 79 | + |
| 80 | + const std::string typeName = CRSLinearSystem::GetClass()->className + CRSLinearSystem::GetCustomTemplateName(); |
| 81 | + |
| 82 | + py::class_<CRSLinearSystem, |
| 83 | + sofa::core::objectmodel::BaseObject, |
| 84 | + sofapython3::py_shared_ptr<CRSLinearSystem> > c(m, typeName.c_str(), sofapython3::doc::linearsystem::linearSystemClass); |
| 85 | + |
| 86 | + c.def("A", [](CRSLinearSystem& self) -> EigenSparseMatrix<Real> |
| 87 | + { |
| 88 | + if (CRS* matrix = self.getSystemMatrix()) |
| 89 | + { |
| 90 | + if (matrix->colsValue.empty()) //null matrix: contains no entries |
| 91 | + { |
| 92 | + return EigenSparseMatrix<Real>{matrix->rows(), matrix->cols()}; |
| 93 | + } |
| 94 | + return toEigen(*matrix); |
| 95 | + } |
| 96 | + return {}; |
| 97 | + }, sofapython3::doc::linearsystem::linearSystem_A); |
| 98 | + |
| 99 | + c.def("b", [](CRSLinearSystem& self) -> Vector<Real> |
| 100 | + { |
| 101 | + if (auto* vector = self.getRHSVector()) |
| 102 | + { |
| 103 | + return EigenVectorMap<Real>(vector->ptr(), vector->size()); |
| 104 | + } |
| 105 | + return {}; |
| 106 | + }, sofapython3::doc::linearsystem::linearSystem_b); |
| 107 | + |
| 108 | + c.def("x", [](CRSLinearSystem& self) -> Vector<Real> |
| 109 | + { |
| 110 | + if (auto* vector = self.getSolutionVector()) |
| 111 | + { |
| 112 | + return EigenVectorMap<Real>(vector->ptr(), vector->size()); |
| 113 | + } |
| 114 | + return {}; |
| 115 | + }, sofapython3::doc::linearsystem::linearSystem_x); |
| 116 | + |
| 117 | + |
| 118 | + /// register the binding in the downcasting subsystem |
| 119 | + PythonFactory::registerType<CRSLinearSystem>([](sofa::core::objectmodel::Base* object) |
| 120 | + { |
| 121 | + return py::cast(dynamic_cast<CRSLinearSystem*>(object)); |
| 122 | + }); |
| 123 | +} |
| 124 | + |
| 125 | +void moduleAddLinearSystem(py::module &m) |
| 126 | +{ |
| 127 | + bindLinearSystems<SReal>(m); |
| 128 | + bindLinearSystems<sofa::type::Mat<3, 3, SReal> >(m); |
| 129 | +} |
| 130 | + |
| 131 | +} |
0 commit comments