imate
C++/CUDA Reference
Loading...
Searching...
No Matches
c_golub_kahn_bidiagonalization.cpp
Go to the documentation of this file.
1/*
2 * SPDX-FileCopyrightText: Copyright 2021, Siavash Ameli <sameli@berkeley.edu>
3 * SPDX-License-Identifier: BSD-3-Clause
4 * SPDX-FileType: SOURCE
5 *
6 * This program is free software: you can redistribute it and/or modify it
7 * under the terms of the license found in the LICENSE.txt file in the root
8 * directory of this source tree.
9 */
10
11
12// =======
13// Headers
14// =======
15
17#include <cmath> // std::sqrt
18#include "./c_orthogonalization.h" // cOrthogonalization
19#include "../_c_basic_algebra/c_vector_operations.h" // cVectorOperations
20
21
22// ============================
23// golub-kahn bidiagonalization
24// ============================
25
109
110template <typename DataType>
113 const DataType* v,
114 const LongIndexType n,
115 const IndexType m,
116 const DataType lanczos_tol,
117 const FlagType orthogonalize,
118 DataType* alpha,
119 DataType* beta)
120{
121 // buffer_size is number of last orthogonal vectors to keep in buffers U, V
122 IndexType buffer_size;
123 if (orthogonalize == 0)
124 {
125 // At least two vectors must be stored in buffer for Lanczos recursion
126 buffer_size = 2;
127 }
128 else if ((orthogonalize < 0) ||
129 (orthogonalize > static_cast<FlagType>(m) - 1))
130 {
131 // Using full reorthogonalization, keep all of the m vectors in buffer
132 buffer_size = m;
133 }
134 else
135 {
136 // Orthogonalize with less than m vectors (0 < orthogonalize < m-1)
137 // plus one vector for the latest (the j-th) vector
138 buffer_size = orthogonalize + 1;
139 }
140
141 // Allocate 2D array (as 1D array, and coalesced row-wise) to store
142 // the last buffer_size of orthogonalized vectors of length n. New vectors
143 // are stored by cycling through the buffer to replace with old ones.
144 DataType* U = new DataType[n * buffer_size];
145 DataType* V = new DataType[n * buffer_size];
146
147 // Normalize vector v and copy to v_old
149
150 // Declare iterators
151 IndexType j;
152 IndexType lanczos_size = 0;
153 IndexType num_ortho;
154
155 // Golub-Kahn iteration
156 for (j=0; j < m; ++j)
157 {
158 // Counter for the non-zero size of alpha and beta
159 ++lanczos_size;
160
161 // u_new = A.dot(v_old)
162 A->dot(&V[(j % buffer_size)*n], &U[(j % buffer_size)*n]);
163
164 // Performing: u_new[i] = u_new[i] - beta[j] * u_old[i]
165 if (j > 0)
166 {
168 &U[((j-1) % buffer_size)*n], n, beta[j-1],
169 &U[(j % buffer_size)*n]);
170 }
171
172 // orthogonalize u_new against previous vectors
173 if (orthogonalize != 0)
174 {
175 // Find how many column vectors are filled so far in the buffer V
176 if (j < buffer_size)
177 {
178 num_ortho = j;
179 }
180 else
181 {
182 num_ortho = buffer_size - 1;
183 }
184
185 // Gram-Schmidt process
186 if (j > 0)
187 {
189 &U[0], n, buffer_size, (j-1)%buffer_size, num_ortho,
190 &U[(j % buffer_size)*n]);
191 }
192 }
193
194 // Normalize u_new and set its norm to alpha[j]
196 &U[(j % buffer_size)*n], n);
197
198 // Performing: v_new = A.T.dot(u_new) - alpha[j] * v_old
199 A->transpose_dot(&U[(j % buffer_size)*n], &V[((j+1) % buffer_size)*n]);
200
201 // Performing: v_new[i] = v_new[i] - alpha[j] * v_old[i]
203 &V[(j % buffer_size)*n], n, alpha[j],
204 &V[((j+1) % buffer_size)*n]);
205
206 // orthogonalize v_new against previous vectors
207 if (orthogonalize != 0)
208 {
210 &V[0], n, buffer_size, j%buffer_size, num_ortho,
211 &V[((j+1) % buffer_size)*n]);
212 }
213
214 // Update beta as the norm of v_new
216 &V[((j+1) % buffer_size)*n], n);
217
218 // Exit criterion when the vector r is zero. If each component of a
219 // zero vector has the tolerance epsilon, (which is called lanczos_tol
220 // here), the tolerance of norm of r is epsilon times sqrt of n.
221 if (beta[j] < lanczos_tol * static_cast<DataType>(std::sqrt(n)))
222 {
223 break;
224 }
225 }
226
227 // Free dynamic memory
228 delete[] U;
229 delete[] V;
230
231 return lanczos_size;
232}
233
234
235// ===============================
236// Explicit template instantiation
237// ===============================
238
239// golub kahn bidiagonalization
242 const float* v,
243 const LongIndexType n,
244 const IndexType m,
245 const float lanczos_tol,
246 const FlagType orthogonalize,
247 float* alpha,
248 float* beta);
249
252 const double* v,
253 const LongIndexType n,
254 const IndexType m,
255 const double lanczos_tol,
256 const FlagType orthogonalize,
257 double* alpha,
258 double* beta);
259
262 const long double* v,
263 const LongIndexType n,
264 const IndexType m,
265 const long double lanczos_tol,
266 const FlagType orthogonalize,
267 long double* alpha,
268 long double* beta);
template IndexType c_golub_kahn_bidiagonalization< float >(cLinearOperator< float > *A, const float *v, const LongIndexType n, const IndexType m, const float lanczos_tol, const FlagType orthogonalize, float *alpha, float *beta)
IndexType c_golub_kahn_bidiagonalization(cLinearOperator< DataType > *A, const DataType *v, const LongIndexType n, const IndexType m, const DataType lanczos_tol, const FlagType orthogonalize, DataType *alpha, DataType *beta)
Bi-diagonalizes the positive-definite matrix A using Golub-Kahn-Lanczos method.
template IndexType c_golub_kahn_bidiagonalization< long double >(cLinearOperator< long double > *A, const long double *v, const LongIndexType n, const IndexType m, const long double lanczos_tol, const FlagType orthogonalize, long double *alpha, long double *beta)
template IndexType c_golub_kahn_bidiagonalization< double >(cLinearOperator< double > *A, const double *v, const LongIndexType n, const IndexType m, const double lanczos_tol, const FlagType orthogonalize, double *alpha, double *beta)
Base class for linear operators. This class serves as interface for all derived classes.
virtual void transpose_dot(const DataType *vector, DataType *product)=0
virtual void dot(const DataType *vector, DataType *product)=0
static void gram_schmidt_process(const DataType *V, const LongIndexType vector_size, const IndexType num_vectors, const IndexType last_vector, const FlagType num_ortho, DataType *r)
Modified Gram-Schmidt orthogonalization process to orthogonalize the vector v against a subset of the...
static DataType normalize_vector_and_copy(const DataType *RESTRICT vector, const LongIndexType vector_size, DataType *RESTRICT output_vector)
Normalizes a vector based on Euclidean 2-norm. The result is written into another vector.
static DataType normalize_vector_in_place(DataType *RESTRICT vector, const LongIndexType vector_size)
Normalizes a vector based on Euclidean 2-norm. The result is written in-place.
static void subtract_scaled_vector(const DataType *RESTRICT input_vector, const LongIndexType vector_size, const DataType scale, DataType *RESTRICT output_vector)
Subtracts the scaled input vector from the output vector.
int LongIndexType
Definition types.h:60
int FlagType
Definition types.h:68
int IndexType
Definition types.h:65