imate
C++/CUDA Reference
c_golub_kahn_bidiagonalization.h File Reference
Include dependency graph for c_golub_kahn_bidiagonalization.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

template<typename DataType >
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. More...
 

Function Documentation

◆ c_golub_kahn_bidiagonalization()

template<typename DataType >
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.

This method bi-diagonalizes matrix A to B using the start vector w. m is the Lanczos degree, which will be the size of square matrix B.

The output of this function are alpha (of length m) and beta (of length m+1) which are diagonal (alpha[:]) and off-diagonal (beta[1:]) elements of the bi-diagonal (m,m) symmetric and positive-definite matrix B.

Lanczos tridiagonalization vs Golub-Kahn Bidiagonalization

  • The Lanczos tri-diagonalization is twice faster (in runtime), as it has only one matrix-vector multiplication. Whereas the Golub-Kahn bi-diagonalization has two matrix-vector multiplications.
  • The Lanczos tri-diagonalization can only be applied to symmetric matrices. Whereas the Golub-Kahn bi-diagonalization can be applied to any matrix.

Reference

  • NetLib Algorithm 6.27, netlib.org/utk/people/JackDongarra/etemplates/node198.html
  • Matrix Computations, Golub, p. 495
  • Demmel, J., Templates for Solution of Algebraic Eigenvalue Problem, p. 143
Warning
When the matrix A is very close to the identity matrix, the Golub-Kahn bi-diagonalization method can not find beta, as beta becomes zero. If A is not exactly identity, you may decrease the Tolerance to a very small number. However, if A is almost identity matrix, decreasing lanczos_tol will not help, and this function cannot be used.
See also
lanczos_tridiagonalizaton
Parameters
[in]AA linear operator that represents a matrix of size (n,n) and can perform matrix-vector operation with dot() method and transposed matrix-vector operation with transpose_dot() method. This matrix should be positive-definite.
[in]vStart vector for the Lanczos tri-diagonalization. Column vector of size n. It could be generated randomly. Often it is generated by the Rademacher distribution with entries +1 and -1.
[in]nSize of the square matrix A, which is also the size of the vector v.
[in]mLanczos degree, which is the number of Lanczos iterations.
[in]lanczos_tolThe tolerance of the residual error of the Lanczos iteration.
[in]orthogonalizeIndicates whether to orthogonalize the orthogonal eigenvectors during Lanczos recursive iterations.
  • If set to 0, no orthogonalization is performed.
  • If set to a negative integer, a newly computed eigenvector is orthogonalized against all the previous eigenvectors (full reorthogonalization).
  • If set to a positive integer, say q less than lanczos_degree, the newly computed eigenvector is orthogonalized against the last q previous eigenvectors (partial reorthogonalization).
  • If set to an integer larger than lanczos_degree, it is cut to lanczos_degree, which effectively orthogonalizes against all previous eigenvectors (full reorthogonalization).
[out]alphaThis is a 1D array of size m and alpha[:] constitute the diagonal elements of the bi-diagonal matrix B. This is the output and written in place.
[out]betaThis is a 1D array of size m, and the elements beta[:] constitute the sup-diagonals of the bi-diagonal matrix B. This array is the output and written in place.
Returns
Counter for the Lanczos iterations. Normally, the size of the output matrix should be (m,m), which is the Lanczos degree. However, if the algorithm terminates early, the size of alpha and beta, and hence the output tri-diagonal matrix, is smaller. This counter keeps track of the non-zero size of alpha and beta.

Definition at line 111 of file c_golub_kahn_bidiagonalization.cpp.

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 * sqrt(n))
222  {
223  break;
224  }
225  }
226 
227  // Free dynamic memory
228  delete[] U;
229  delete[] V;
230 
231  return lanczos_size;
232 }
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 void subtract_scaled_vector(const DataType *input_vector, const LongIndexType vector_size, const DataType scale, DataType *output_vector)
Subtracts the scaled input vector from the output vector.
static DataType normalize_vector_in_place(DataType *vector, const LongIndexType vector_size)
Normalizes a vector based on Euclidean 2-norm. The result is written in-place.
static DataType normalize_vector_and_copy(const DataType *vector, const LongIndexType vector_size, DataType *output_vector)
Normalizes a vector based on Euclidean 2-norm. The result is written into another vector.
int FlagType
Definition: types.h:68
int IndexType
Definition: types.h:65

References cLinearOperator< DataType >::dot(), cOrthogonalization< DataType >::gram_schmidt_process(), cVectorOperations< DataType >::normalize_vector_and_copy(), cVectorOperations< DataType >::normalize_vector_in_place(), cVectorOperations< DataType >::subtract_scaled_vector(), and cLinearOperator< DataType >::transpose_dot().

Referenced by cTraceEstimator< DataType >::_c_stochastic_lanczos_quadrature().

Here is the call graph for this function:
Here is the caller graph for this function: