forked from yixuan/spectra
-
Notifications
You must be signed in to change notification settings - Fork 0
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
add Gram-Schmidt process #11 #15
Merged
Merged
Changes from 15 commits
Commits
Show all changes
21 commits
Select commit
Hold shift + click to select a range
6b22c19
add Gram-Schmidt process #11
felipeZ 1a4d669
add partial orthonormalization test
felipeZ 730695d
improved gram-schmidt #11
felipeZ b8491e3
fixed the GramSchmidt process #11
felipeZ 761d14f
changed check orthonormal signature #11
felipeZ 043ac01
renamed gramschmidt to orthogonalization
felipeZ 0fef8ba
implemented QR and modified GS #11
felipeZ 06967e9
Merge branch 'master' into gramschmidt
felipeZ dd636d9
renamed twice gs #11
felipeZ 86cd4d1
implement modified Gram-Schmidt #11
felipeZ 539cfc3
used vectorize modified GS #11
felipeZ a7f645f
turned class into free functions
JensWehner 5222e86
rewrote ortho and introduced JW
NicoRenaud 73cae07
fixed JW
NicoRenaud 7249399
removed unnecessary copy, fixed spelling
JensWehner dda0164
Merge branch 'master' into gramschmidt
JensWehner f166d74
make format
JensWehner 45dbd6d
Merge branch 'gramschmidt' of https://github.com/NLESC-JCER/spectra i…
JensWehner 1a9ded2
split check
NicoRenaud 17b8ea1
Merge branch 'gramschmidt' of https://github.com/NLESC-JCER/Spectra i…
NicoRenaud 3881b68
fix small pr from jw
NicoRenaud File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -4,3 +4,5 @@ include/Eigen/* | |
.idea/ | ||
cmake-build*/ | ||
build/* | ||
.vscode | ||
*~ |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,98 @@ | ||
// Copyright (C) 2020 Netherlands eScience Center <[email protected]> | ||
// | ||
// This Source Code Form is subject to the terms of the Mozilla | ||
// Public License v. 2.0. If a copy of the MPL was not distributed | ||
// with this file, You can obtain one at https://mozilla.org/MPL/2.0/. | ||
|
||
#ifndef SPECTRA_ORTHOGONALIZATION_H | ||
#define SPECTRA_ORTHOGONALIZATION_H | ||
|
||
#include <Eigen/Core> | ||
#include <Eigen/Dense> | ||
|
||
namespace Spectra { | ||
|
||
template <typename Matrix> | ||
Eigen::Index sanity_check(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
{ | ||
assert(in_output.cols() > leftColsToSkip && "leftColsToSkip is larger than columns of matrix"); | ||
assert(leftColsToSkip >= 0 && "leftColsToSkip is negative"); | ||
if (leftColsToSkip == 0) | ||
{ | ||
in_output.col(0).normalize(); | ||
leftColsToSkip = 1; | ||
} | ||
return leftColsToSkip; | ||
} | ||
NicoRenaud marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
template <typename Matrix> | ||
void QR_orthogonalisation(Matrix& in_output) | ||
{ | ||
using InternalMatrix= Eigen::Matrix<typename Matrix::Scalar,Eigen::Dynamic, Eigen::Dynamic>; | ||
Eigen::Index nrows = in_output.rows(); | ||
Eigen::Index ncols = in_output.cols(); | ||
ncols = std::min(nrows, ncols); | ||
InternalMatrix I = InternalMatrix::Identity(nrows, ncols); | ||
Eigen::HouseholderQR<Matrix> qr(in_output); | ||
in_output = qr.householderQ() * I; | ||
} | ||
|
||
template <typename Matrix> | ||
void MGS_orthogonalisation(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
{ | ||
leftColsToSkip = sanity_check(in_output, leftColsToSkip); | ||
|
||
for (Eigen::Index k = leftColsToSkip; k < in_output.cols(); ++k) | ||
{ | ||
for (Eigen::Index j = 0; j < k; j++) | ||
{ | ||
in_output.col(k) -= in_output.col(j).dot(in_output.col(k)) / (in_output.col(j).dot(in_output.col(j))) * in_output.col(j); | ||
NicoRenaud marked this conversation as resolved.
Show resolved
Hide resolved
|
||
} | ||
in_output.col(k).normalize(); | ||
} | ||
} | ||
|
||
template <typename Matrix> | ||
void GS_orthogonalisation(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
{ | ||
leftColsToSkip = sanity_check(in_output, leftColsToSkip); | ||
|
||
for (Eigen::Index j = leftColsToSkip; j < in_output.cols(); ++j) | ||
{ | ||
in_output.col(j) -= in_output.leftCols(j) * (in_output.leftCols(j).transpose() * in_output.col(j)); | ||
in_output.col(j).normalize(); | ||
} | ||
} | ||
|
||
template <typename Matrix> | ||
void twice_is_enough_orthogonalisation(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
{ | ||
GS_orthogonalisation(in_output, leftColsToSkip); | ||
GS_orthogonalisation(in_output, leftColsToSkip); | ||
} | ||
|
||
template <typename Matrix> | ||
void partial_orthogonalisation(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
NicoRenaud marked this conversation as resolved.
Show resolved
Hide resolved
|
||
{ | ||
leftColsToSkip = sanity_check(in_output, leftColsToSkip); | ||
|
||
Eigen::Index rightColToOrtho = in_output.cols() - leftColsToSkip; | ||
in_output.rightCols(rightColToOrtho) -= (in_output.leftCols(leftColsToSkip) * in_output.leftCols(leftColsToSkip).transpose()) * in_output.rightCols(rightColToOrtho); | ||
in_output.rightCols(rightColToOrtho).colwise().normalize(); | ||
} | ||
|
||
template <typename Matrix> | ||
void JensWehner_orthogonalisation(Matrix& in_output, Eigen::Index leftColsToSkip = 0) | ||
{ | ||
leftColsToSkip = sanity_check(in_output, leftColsToSkip); | ||
|
||
Eigen::Index rightColToOrtho = in_output.cols() - leftColsToSkip; | ||
partial_orthogonalisation(in_output, leftColsToSkip); | ||
Eigen::Ref<Matrix> right_cols = in_output.rightCols(rightColToOrtho); | ||
QR_orthogonalisation(right_cols); | ||
in_output.rightCols(rightColToOrtho) = right_cols; | ||
} | ||
|
||
} // namespace Spectra | ||
|
||
#endif //SPECTRA_ORTHOGONALIZATION_H |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,99 @@ | ||
#include <Eigen/Core> | ||
#include <Spectra/LinAlg/Orthogonalization.h> | ||
#include <iostream> | ||
using namespace Spectra; | ||
|
||
#define CATCH_CONFIG_MAIN | ||
#include "catch.hpp" | ||
|
||
using Eigen::MatrixXd; | ||
using Eigen::VectorXd; | ||
using Eigen::Index; | ||
|
||
template <typename Matrix> | ||
void check_orthogonality(const Matrix& basis) | ||
{ | ||
const double tol = 1e-12; | ||
Matrix xs = basis.transpose() * basis; | ||
INFO("The orthonormalized basis must fulfill that basis.T * basis = I"); | ||
INFO("Matrix is\n " << basis); | ||
INFO("Overlap is\n " << xs); | ||
CHECK(xs.isIdentity(tol)); | ||
} | ||
|
||
TEST_CASE("complete orthonormalization", "[orthogonalisation]") | ||
{ | ||
std::srand(123); | ||
const Index n = 20; | ||
|
||
MatrixXd mat = MatrixXd::Random(n, n); | ||
|
||
SECTION("MGS") | ||
{ | ||
MGS_orthogonalisation(mat); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("GS") | ||
{ | ||
GS_orthogonalisation(mat); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("QR") | ||
{ | ||
QR_orthogonalisation(mat); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("twice_is_enough") | ||
{ | ||
twice_is_enough_orthogonalisation(mat); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("JensWehner") | ||
{ | ||
JensWehner_orthogonalisation(mat); | ||
check_orthogonality(mat); | ||
} | ||
} | ||
|
||
TEST_CASE("Partial orthonormalization", "[orthogonalisation]") | ||
{ | ||
std::srand(123); | ||
const Index n = 20; | ||
const Index sub = 5; | ||
Index start = n - sub; | ||
|
||
// Create a n x 20 orthonormal basis | ||
MatrixXd mat = MatrixXd::Random(n, start); | ||
QR_orthogonalisation(mat); | ||
|
||
mat.conservativeResize(Eigen::NoChange, n); | ||
mat.rightCols(sub) = MatrixXd::Random(n, sub); | ||
|
||
SECTION("MGS") | ||
{ | ||
MGS_orthogonalisation(mat, start); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("GS") | ||
{ | ||
GS_orthogonalisation(mat, start); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("twice_is_enough") | ||
{ | ||
twice_is_enough_orthogonalisation(mat, start); | ||
check_orthogonality(mat); | ||
} | ||
|
||
SECTION("JensWehner") | ||
{ | ||
JensWehner_orthogonalisation(mat, start); | ||
check_orthogonality(mat); | ||
} | ||
} |
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
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.
#include <Eigen/QR>
seems sufficient.#include <Eigen/Dense>
will also include many heavy solvers such as LU, SVD, etc.