-
Notifications
You must be signed in to change notification settings - Fork 65
adding weights for C++ (#267) #268
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
base: main
Are you sure you want to change the base?
Changes from all commits
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -18,8 +18,9 @@ | |
| #include "utils.h" | ||
| #include <cassert> | ||
|
|
||
| #ifdef EIGEN | ||
| //#ifdef EIGEN | ||
| #include <eigen3/Eigen/SparseLU> | ||
| #include <eigen3/Eigen/SparseQR> | ||
|
|
||
| vec Utils::spsolve_eigen(const sp_mat &A, const vec &b) { | ||
| Eigen::SparseMatrix<Real> eigen_A(A.n_rows, A.n_cols); | ||
|
|
@@ -47,7 +48,34 @@ vec Utils::spsolve_eigen(const sp_mat &A, const vec &b) { | |
|
|
||
| return vec(eigen_x.data(), eigen_x.size()); | ||
| } | ||
| #endif | ||
|
|
||
| vec Utils::spsolve_eigenQR(const sp_mat &A, const vec &b) { | ||
|
Collaborator
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Could you elaborate on why QR factorization is needed here as opposed to LU? I’m assuming this is for computing P and Q. Are we expecting any ill-conditioned matrices?
Collaborator
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Moreover, the only difference between the LU method and the QR one is: 1c1
< vec Utils::spsolve_eigen(const sp_mat &A, const vec &b) {
---
> vec Utils::spsolve_eigenQR(const sp_mat &A, const vec &b) {
4c4
< Eigen::SparseLU<Eigen::SparseMatrix<Real>, Eigen::COLAMDOrdering<int>> solver;
---
> Eigen::SparseQR<Eigen::SparseMatrix<Real>, Eigen::COLAMDOrdering<int>> solver;We should consider refactoring to avoid code duplication. |
||
| Eigen::SparseMatrix<Real> eigen_A(A.n_rows, A.n_cols); | ||
| std::vector<Eigen::Triplet<Real>> triplets; | ||
| Eigen::SparseQR<Eigen::SparseMatrix<Real>, Eigen::COLAMDOrdering<int>> solver; | ||
|
|
||
| Eigen::VectorXd eigen_x(A.n_rows); | ||
| triplets.reserve(5 * A.n_rows); | ||
|
|
||
| auto it = A.begin(); | ||
| while (it != A.end()) { | ||
| triplets.push_back(Eigen::Triplet<Real>(it.row(), it.col(), *it)); | ||
| ++it; | ||
| } | ||
|
|
||
| eigen_A.setFromTriplets(triplets.begin(), triplets.end()); | ||
| triplets.clear(); | ||
|
|
||
| auto b_ = conv_to<std::vector<Real>>::from(b); | ||
| Eigen::Map<Eigen::VectorXd> eigen_b(b_.data(), b_.size()); | ||
|
|
||
| solver.analyzePattern(eigen_A); | ||
| solver.factorize(eigen_A); | ||
| eigen_x = solver.solve(eigen_b); | ||
|
|
||
| return vec(eigen_x.data(), eigen_x.size()); | ||
| } | ||
| //#endif | ||
|
|
||
| // Basic implementation of Kronecker product | ||
| /* | ||
|
|
||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,32 @@ | ||
| /* | ||
| * SPDX-License-Identifier: GPL-3.0-or-later | ||
| * © 2008-2024 San Diego State University Research Foundation (SDSURF). | ||
| * See LICENSE file or https://www.gnu.org/licenses/gpl-3.0.html for details. | ||
| */ | ||
|
|
||
| /* | ||
| * @file weights.cpp | ||
| * | ||
| * @brief Mimetic Differences Q Weights | ||
| * | ||
| * @date 2025/10/14 | ||
| * | ||
| */ | ||
|
|
||
| #include "utils.h" | ||
| #include "divergence.h" | ||
| #include "gradient.h" | ||
| #include "weightsP.h" | ||
|
|
||
| WeightsP::WeightsP(u16 k, u32 m, Real dx) | ||
| { | ||
| Gradient G(k, m, dx); | ||
|
|
||
| arma::vec b(m+2); | ||
| b.at(0) = -1.0; | ||
| b.at(m+1) = 1.0; | ||
|
|
||
|
Collaborator
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. nit: extra lineIDK |
||
| sp_mat Gtranspose = G.t(); | ||
|
Collaborator
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. arma::sp_mat |
||
| *this = Utils::spsolve_eigenQR(Gtranspose,b); | ||
|
|
||
| } | ||
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,42 @@ | ||
| /* | ||
| * SPDX-License-Identifier: GPL-3.0-or-later | ||
| * © 2008-2024 San Diego State University Research Foundation (SDSURF). | ||
| * See LICENSE file or https://www.gnu.org/licenses/gpl-3.0.html for details. | ||
| */ | ||
|
|
||
| /* | ||
| * @file weights.h | ||
| * | ||
| * @brief Mimetic Differences Weights | ||
| * | ||
| * @date 2025/10/14 | ||
| * | ||
| */ | ||
|
|
||
| #ifndef WEIGHTSP_H | ||
| #define WEIGHTSP_H | ||
|
|
||
| #include "utils.h" | ||
| #include <cassert> | ||
|
|
||
| /** | ||
| * @brief Generate P Weights | ||
| * | ||
| */ | ||
| class WeightsP : public sp_vec { | ||
|
|
||
| public: | ||
| using sp_vec::operator=; | ||
|
|
||
| /** | ||
| * @brief P Weights Constructor | ||
| * | ||
| * @param k Order of accuracy | ||
| * @param m Number of cells | ||
| * @param dx Spacing between cells | ||
| */ | ||
| WeightsP(u16 k, u32 m, Real dx); | ||
|
|
||
| }; | ||
|
|
||
| #endif // WEIGHTSP_H |
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,43 @@ | ||
| /* | ||
| * SPDX-License-Identifier: GPL-3.0-or-later | ||
| * © 2008-2024 San Diego State University Research Foundation (SDSURF). | ||
| * See LICENSE file or https://www.gnu.org/licenses/gpl-3.0.html for details. | ||
| */ | ||
|
|
||
| /* | ||
| * @file weights.cpp | ||
| * | ||
| * @brief Mimetic Differences Q Weights | ||
| * | ||
| * @date 2025/10/14 | ||
| * | ||
| */ | ||
|
|
||
| #include "utils.h" | ||
| #include "divergence.h" | ||
| #include "gradient.h" | ||
| #include "weightsQ.h" | ||
|
|
||
| WeightsQ::WeightsQ(u16 k, u32 m, Real dx): sp_vec(m + 1) | ||
| { | ||
| Divergence D(k, m, dx); | ||
| D.shed_row(0); | ||
| D.shed_row(m); | ||
|
|
||
| vec b(m+1); | ||
| b.at(0) = -1.0; | ||
| b.at(m) = 1.0; | ||
|
|
||
| sp_mat Dtranspose = D.t(); | ||
| sp_vec Q_prime = Utils::spsolve_eigenQR(Dtranspose, b); | ||
| sp_vec Q(m+2); | ||
| Q.at(0) = 1.0; | ||
| Q.at(m+1) = 1.0; | ||
| for (int i = 1; i < m+1; ++i) { | ||
| Q.at(i) = Q_prime.at(i-1); | ||
| } | ||
|
|
||
| *this = Q; | ||
|
|
||
| } | ||
|
|
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,42 @@ | ||
| /* | ||
| * SPDX-License-Identifier: GPL-3.0-or-later | ||
| * © 2008-2024 San Diego State University Research Foundation (SDSURF). | ||
| * See LICENSE file or https://www.gnu.org/licenses/gpl-3.0.html for details. | ||
| */ | ||
|
|
||
| /* | ||
| * @file weights.h | ||
| * | ||
| * @brief Mimetic Differences Weights | ||
| * | ||
| * @date 2025/10/14 | ||
| * | ||
| */ | ||
|
|
||
| #ifndef WEIGHTSQ_H | ||
| #define WEIGHTSQ_H | ||
|
|
||
| #include "utils.h" | ||
| #include <cassert> | ||
|
|
||
| /** | ||
| * @brief Generate Q Weights | ||
| * | ||
| */ | ||
| class WeightsQ : public sp_vec { | ||
|
|
||
| public: | ||
| using sp_vec::operator=; | ||
|
|
||
| /** | ||
| * @brief Q Weights Constructor | ||
| * | ||
| * @param k Order of accuracy | ||
| * @param m Number of cells | ||
| * @param dx Spacing between cells | ||
| */ | ||
| WeightsQ(u16 k, u32 m, Real dx); | ||
|
|
||
| }; | ||
|
|
||
| #endif // WEIGHTSQ_H |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Why are we commenting out this line? Is Eigen mandatory now?