mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
gramSchmidt.hpp
Go to the documentation of this file.
1/** \file gramSchmidt.hpp
2 * \brief Procedures to orthogonalize vector basis sets
3 *
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 * \ingroup signal_processing_files
7 *
8 */
9
10//***********************************************************************//
11// Copyright 2015, 2016, 2017 Jared R. Males (jaredmales@gmail.com)
12//
13// This file is part of mxlib.
14//
15// mxlib is free software: you can redistribute it and/or modify
16// it under the terms of the GNU General Public License as published by
17// the Free Software Foundation, either version 3 of the License, or
18// (at your option) any later version.
19//
20// mxlib is distributed in the hope that it will be useful,
21// but WITHOUT ANY WARRANTY; without even the implied warranty of
22// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23// GNU General Public License for more details.
24//
25// You should have received a copy of the GNU General Public License
26// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
27//***********************************************************************//
28
29#ifndef gramSchmidt_hpp
30#define gramSchmidt_hpp
31
32namespace mx
33{
34namespace sigproc
35{
36
37/// Perform Gram-Schmidt ortogonalization of a basis set, and normalize the result.
38/** Performs the stabilized Gram-Schmidt procedure on the input basis set, which
39 * is the columns of the array. Optionally the output is normalized.
40 *
41 * \tparam progress if true, then the loop index is printed for progress reporting
42 * \tparam eigenTout is the Eigen array type of the desired output
43 * \tparam eigenTin is the Eigen array type of the input
44 *
45 * \ingroup signal_processing
46 */
47template <int progress = 0, typename eigenTout, typename eigenTin>
48void gramSchmidt( eigenTout &out, ///< [out] the orthonormal basis set constructed from the input
49 const eigenTin &in, ///< [in] a basis set, where each column represents one vector.
50 bool normalize = true ///< [in] [optional] whether or not to normalize the output
51)
52{
53 out.resize( in.rows(), in.cols() );
54
55 out.col( 0 ) = in.col( 0 );
56
57 for( int i = 1; i < in.cols(); ++i )
58 {
59 if( progress )
60 {
61 std::cout << i + 1 << "/" << in.cols() << "\n";
62 }
63
64 // out.col(i) = in.col(i);
65
66 out.col( i ) = in.col( i ) - ( ( in.col( i ).matrix().dot( out.col( 0 ).matrix() ) ) /
67 ( out.col( 0 ).matrix().dot( out.col( 0 ).matrix() ) ) ) *
68 out.col( 0 );
69
70 for( int j = 1; j < i; ++j )
71 {
72 out.col( i ) = out.col( i ) - ( ( out.col( i ).matrix().dot( out.col( j ).matrix() ) ) /
73 ( out.col( j ).matrix().dot( out.col( j ).matrix() ) ) ) *
74 out.col( j );
75 }
76 }
77
78 if( normalize )
79 {
80 for( int i = 0; i < out.cols(); ++i )
81 {
82 out.col( i ) = out.col( i ) / out.col( i ).matrix().norm();
83 }
84 }
85}
86
87/// Perform Gram-Schmidt ortogonalization of a basis set on a window, and normalize the result.
88/** Performs the stabilized Gram-Schmidt procedure on the input basis set over a window (or
89 * weight function), followed by normalization of the result.
90 *
91 * \param out [out] is the orthonormal basis set constructed from the input
92 * \param in [in] is a basis set, where each column represents one vector.
93 * \param window [in] is the window, or weighting function
94 *
95 * \tparam progress if true, then the loop index is printed for progress reporting
96 * \tparam eigenTout is the Eigen array type of the desired output
97 * \tparam eigenTin is the Eigen array type of the input
98 * \tparam eigenTWin is the Eigen array type of the window
99 *
100 * \ingroup signal_processing
101 */
102template <int progress = 0, typename eigenTout, typename eigenTin, typename eigenTWin>
103void gramSchmidt( eigenTout &out, const eigenTin &in, const eigenTWin &window )
104{
105 // out.resize(in.rows(), in.cols());
106
107 out.col( 0 ) = in.col( 0 );
108
109 for( int i = 1; i < in.cols(); ++i )
110 {
111 if( progress )
112 {
113 std::cout << i + 1 << "/" << in.cols() << "\n";
114 }
115
116 // out.col(i) = in.col(i);
117
118 out.col( i ) = in.col( i ) - ( ( ( in.col( i ) * window ).matrix().dot( out.col( 0 ).matrix() ) ) /
119 ( ( out.col( 0 ) * window ).matrix().dot( out.col( 0 ).matrix() ) ) ) *
120 out.col( 0 );
121
122 for( int j = 1; j < i; ++j )
123 {
124 out.col( i ) = out.col( i ) - ( ( ( out.col( i ) * window ).matrix().dot( out.col( j ).matrix() ) ) /
125 ( ( out.col( j ) * window ).matrix().dot( out.col( j ).matrix() ) ) ) *
126 out.col( j );
127 }
128 }
129
130 for( int i = 0; i < out.cols(); ++i )
131 {
132 out.col( i ) = out.col( i ) / ( out.col( i ) * window.sqrt() ).matrix().norm();
133 }
134}
135
136// Unwraps the gram schmidt coefficients to give a spectrum in terms of the original basis set
137// Helper function for gramSchmidtSpectrum (below)
138template <int progress = 0, typename eigenT>
139void baseSpectrum( eigenT &bspect, eigenT &gsspect )
140{
141 bspect.resize( gsspect.rows(), gsspect.cols() );
142 bspect.setZero();
143
144 // #pragma omp parallel for
145 for( int i = 0; i < gsspect.rows(); ++i )
146 {
147 bspect( i, i ) = gsspect( i, i );
148
149 for( int j = i - 1; j >= 0; --j )
150 {
151 bspect.row( i ) -= gsspect( i, j ) * bspect.row( j );
152 }
153 }
154}
155
156/// Perform Gram-Schmidt ortogonalization of a basis set, and normalize the result, while recording the spectrum.
157/** Performs the stabilized Gram-Schmidt procedure on the input basis set, followed
158 * by normalization of the result. Also records the spectrum, that is the coefficients of the linear expansion
159 * in the orginal basis set for the resultant basis set.
160 *
161 *
162 * \tparam progress if true, then the loop index is printed for progress reporting
163 * \tparam eigenTout is the Eigen array type of the output orthogonalized array
164 * \tparam eigenTout2 is the Eigen array type of the spectrum
165 * \tparam eigenTin is the Eigen array type of the input
166 *
167 * \ingroup signal_processing
168 */
169template <int progress = 0, typename eigenTout, typename eigenTout2, typename eigenTin>
171 eigenTout &out, ///< [out] the orthonormal basis set constructed from the input
172 eigenTout2 &spect, ///< [out] the spectrum
173 const eigenTin &in, ///< [in] a basis set, where each column represents one vector
174 typename eigenTin::Scalar normPix = 0.0 /**< [in] [optional] area of (usually number of pixels in) the orthogonal
175 region for normalization. If 0 the basis is not renormalized */
176)
177{
178 typedef typename eigenTout::Scalar Scalar;
179
180 out.resize( in.rows(), in.cols() );
181
182 eigenTout2 gsspect;
183
184 gsspect.resize( in.cols(), in.cols() );
185 gsspect.setZero();
186
187 out.col( 0 ) = in.col( 0 );
188 gsspect( 0, 0 ) = 1;
189
190 for( int i = 1; i < in.cols(); ++i )
191 {
192 if( progress )
193 {
194 std::cout << i + 1 << "/" << in.cols() << "\n";
195 }
196
197 gsspect( i, i ) = 1;
198
199 gsspect( i, 0 ) = ( ( in.col( i ).matrix().dot( out.col( 0 ).matrix() ) ) /
200 ( out.col( 0 ).matrix().dot( out.col( 0 ).matrix() ) ) );
201 out.col( i ) = in.col( i ) - gsspect( i, 0 ) * out.col( 0 );
202
203 for( int j = 1; j < i; ++j )
204 {
205 gsspect( i, j ) = ( ( out.col( i ).matrix().dot( out.col( j ).matrix() ) ) /
206 ( out.col( j ).matrix().dot( out.col( j ).matrix() ) ) );
207 out.col( i ) = out.col( i ) - gsspect( i, j ) * out.col( j );
208 }
209 }
210
211 // Here we unwrap the gram schmidt coefficients, giving us the coefficients in the original basis
212 baseSpectrum( spect, gsspect );
213
214 if( normPix > 0 )
215 {
216 Scalar norm;
217
218 for( int i = 0; i < out.cols(); ++i )
219 {
220 norm = sqrt( out.col( i ).square().sum() / normPix );
221
222 out.col( i ) /= norm;
223 // spect.row(i) /= norm;
224 }
225 }
226}
227
228} // namespace sigproc
229} // namespace mx
230
231#endif // gramSchmidt_hpp
void gramSchmidtSpectrum(eigenTout &out, eigenTout2 &spect, const eigenTin &in, typename eigenTin::Scalar normPix=0.0)
Perform Gram-Schmidt ortogonalization of a basis set, and normalize the result, while recording the s...
void gramSchmidt(eigenTout &out, const eigenTin &in, bool normalize=true)
Perform Gram-Schmidt ortogonalization of a basis set, and normalize the result.
The mxlib c++ namespace.
Definition mxError.hpp:106