imate
C++/CUDA Reference
cOrthogonalization< DataType > Class Template Reference

A static class for orthogonalization of vector bases. This class acts as a templated namespace, where all member methods are public and static. More...

#include <c_orthogonalization.h>

Static Public Member Functions

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 column vectors in the array V. More...
 
static void orthogonalize_vectors (DataType *vectors, const LongIndexType vector_size, const IndexType num_vectors)
 Orthogonalizes set of vectors mutually using modified Gram-Schmidt process. More...
 

Detailed Description

template<typename DataType>
class cOrthogonalization< DataType >

A static class for orthogonalization of vector bases. This class acts as a templated namespace, where all member methods are public and static.

See also
RandomVectors

Definition at line 35 of file c_orthogonalization.h.

Member Function Documentation

◆ gram_schmidt_process()

template<typename DataType >
void cOrthogonalization< DataType >::gram_schmidt_process ( const DataType *  V,
const LongIndexType  vector_size,
const IndexType  num_vectors,
const IndexType  last_vector,
const FlagType  num_ortho,
DataType *  v 
)
static

Modified Gram-Schmidt orthogonalization process to orthogonalize the vector v against a subset of the column vectors in the array V.

V is 1D array of the length vector_size*num_vectors to represent a 2D array of a set of num_vectors column vectors, each of the length vector_size. The length of v is also vector_size.

v is orthogonalized against the last num_ortho columns of V starting from the column vector of the index last_vector. If the backward indexing from last_vector becomes a negative index, the index wraps around from the last column vector index, i.e., num_vectors-1 .

  • If num_ortho is zero, or if num_vectors is zero, no orthogonalization is performed.
  • If num_ortho is negative (usually set to -1), then v is orthogonalized against all column vectors of V.
  • If num_ortho is larger than num_vectors, then v is orthogonalized against all column vectors of V.
  • If num_ortho is smaller than num_vectors, then v is orthogonalized against the last num_ortho column vectors of V, starting from the column vector with the index last_vector toward its previous vectors. If the iteration runs into negativen column indices, the column indexing wraps around from the end of the columns from num_vectors-1.

The result of the newer v is written in-place in v.

If vector v is identical to one of the vectors in V, the orthogonalization against the identical vector is skipped.

If one of the column vectors of V is zero (have zero norm), that vector is ignored.

Note
It is assumed that the caller function fills the column vectors of V periodically in a wrapped around order from column index 0,1,... to num_vectors-1, and newer vectors are replaced on the wrapped index starting from index 0,1,... again. Thus, V only stores the last num_vectors column vectors. The index of the last filled vector is indicated by last_vector.
Warning
The vector v can be indeed one of the columns of V itself. However, in this case, vector v must NOT be orthogonalized against itself, rather, it should only be orthogonalized to the other vectors in V. For instance, if num_vectors=10, and v is the 3rd vector of V, and if num_ortho is 6, then we may set last_vector=2. Then v is orthogonalized againts the six columns 2,1,0,9,8,7, where the last three of them are wrapped around the end of the columns.
See also
c_golub_kahn_bidiagonalizaton, c_lanczos_bidiagonalization
Parameters
[in]V1D coalesced array of vectors representing a 2D array. The length of this 1D array is vector_size*num_vectors, which indicates a 2D array with the shape (vector_size,num_vectors).
[in]vector_sizeThe length of each vector. If we assume V indicates a 2D vector, this is the number of rows of V.
[in]num_vectorsThe number of column vectors. If we assume V indicates a 2D vector, this the number of columns of V.
[in]last_vectorThe column vectors of the array V are rewritten by the caller function in wrapped-around order. That is, once all the columns (from the zeroth to the num_vector-1 vector) are filled, the next vector is rewritten in the place of the zeroth vector, and the indices of newer vectors wrap around the columns of V. Thus, V only retains the last num_vectors vectors. The column index of the last written vector is given by last_vector. This index is a number between 0 and num_vectors-1. The index of the last i-th vector is winding back from the last vector by last_vector-i+1 mod num_vectors.
[in]num_orthoThe number of vectors to be orthogonalized starting from the last vector. 0 indicates no orthogonalization will be performed and the function just returns. A negative value means all vectors will be orthogonalized. A poisitive value will orthogonalize the given number of vectors. This value cannot be larger than the number of vectors.
[in,out]vThe vector that will be orthogonalized against the columns of V. The length of v is vector_size. This vector is modified in-place.

Definition at line 125 of file c_orthogonalization.cpp.

132 {
133  // Determine how many previous vectors to orthogonalize against
134  IndexType num_steps;
135  if ((num_ortho == 0) || (num_vectors < 2))
136  {
137  // No orthogonalization is performed
138  return;
139  }
140  else if ((num_ortho < 0) ||
141  (num_ortho > static_cast<FlagType>(num_vectors)))
142  {
143  // Orthogonalize against all vectors
144  num_steps = num_vectors;
145  }
146  else
147  {
148  // Orthogonalize against only the last num_ortho vectors
149  num_steps = num_ortho;
150  }
151 
152  // Vectors can be orthogonalized at most to the full basis of the vector
153  // space. Thus, num_steps cannot be larger than the dimension of vector
154  // space, which is vector_size.
155  if (num_steps > static_cast<IndexType>(vector_size))
156  {
157  num_steps = vector_size;
158  }
159 
160  IndexType i;
161  DataType inner_prod;
162  DataType norm;
163  DataType norm_v;
164  DataType epsilon = std::numeric_limits<DataType>::epsilon();
165  DataType distance2;
166 
167  // Iterate over vectors
168  for (IndexType step=0; step < num_steps; ++step)
169  {
170  // i is the index of a column vector in V to orthogonalize v against it
171  if ((last_vector % num_vectors) >= step)
172  {
173  i = (last_vector % num_vectors) - step;
174  }
175  else
176  {
177  // Wrap around negative indices from the end of column index
178  i = (last_vector % num_vectors) - step + num_vectors;
179  }
180 
181  // Norm of j-th vector
183  &V[vector_size*i], vector_size);
184 
185  // Check norm
186  if (norm < epsilon * sqrt(vector_size))
187  {
188  std::cerr << "WARNING: norm of the given vector is too small. " \
189  << "Cannot orthogonalize against zero vector. " \
190  << "Skipping." << std::endl;
191  continue;
192  }
193 
194  // Projection
196  &V[vector_size*i], v, vector_size);
197 
198  // scale for subtraction
199  DataType scale = inner_prod / (norm * norm);
200 
201  // If scale is 1, it is possible that vector v and j-th vector are
202  // identical (or close).
203  if (std::abs(std::abs(scale) - 1.0) <= 2.0 * epsilon)
204  {
205  // Norm of the vector v
207  v, vector_size);
208 
209  // Compute distance between the j-th vector and vector v
210  distance2 = norm_v*norm_v - 2.0*inner_prod + norm*norm;
211 
212  // If distance is zero, do not reorthogonalize i-th against
213  // the j-th vector.
214  if (distance2 < 2.0 * epsilon * vector_size)
215  {
216  continue;
217  }
218  }
219 
220  // Subtraction
222  &V[vector_size*i], vector_size, scale, v);
223  }
224 }
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 inner_product(const DataType *vector1, const DataType *vector2, const LongIndexType vector_size)
Computes Euclidean inner product of two vectors.
static DataType euclidean_norm(const DataType *vector, const LongIndexType vector_size)
Computes the Euclidean norm of a 1D array.
int FlagType
Definition: types.h:68
int IndexType
Definition: types.h:65

References cVectorOperations< DataType >::euclidean_norm(), cVectorOperations< DataType >::inner_product(), and cVectorOperations< DataType >::subtract_scaled_vector().

Referenced by c_golub_kahn_bidiagonalization(), and c_lanczos_tridiagonalization().

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

◆ orthogonalize_vectors()

template<typename DataType >
void cOrthogonalization< DataType >::orthogonalize_vectors ( DataType *  vectors,
const LongIndexType  vector_size,
const IndexType  num_vectors 
)
static

Orthogonalizes set of vectors mutually using modified Gram-Schmidt process.

Note
Let m be the number of vectors (num_vectors), and let n be the size of each vector (vector_size). In general, n is much larger (large matrix size), and m is small, in order of a couple of hundred. But for small matrices (where n could be smaller then m), then each vector can be orthogonalized at most to n other vectors. This is because the dimension of the vector space is n. Thus, if there are extra vectors, each vector is orthogonalized to window of the previous n vector.

If one of the column vectors is identical to one of other column vectors in V, one of the vectors is regenerated by random array and the orthogonalization is repeated.

Note
If two vectors are identical (or the norm of their difference is very small), they cannot be orthogonalized against each other. In this case, one of the vectors is re-generated by new random numbers.
Warning
if num_vectors is larger than vector_size, the orthogonalization fails since not all vectors are independent, and at least one vector becomes zero.
Parameters
[in,out]vectors2D array of size vector_size*num_vectors. This array will be modified in-place and will be output of this function. Note that this is Fortran ordering, meaning that the first index is contiguous. Hence, to call the j-th element of the i-th vector, use &vectors[i*vector_size + j].
[in]num_vectorsNumber of columns of vectors array.
[in]vector_sizeNumber of rows of vectors array.

Definition at line 269 of file c_orthogonalization.cpp.

273 {
274  // Do nothing if there is only one vector
275  if (num_vectors < 2)
276  {
277  return;
278  }
279 
280  IndexType i = 0;
281  IndexType j;
282  IndexType start_j;
283  DataType inner_prod;
284  DataType norm_j;
285  DataType norm_i;
286  DataType epsilon = std::numeric_limits<DataType>::epsilon();
287  IndexType success = 1;
288  IndexType max_num_trials = 20;
289  IndexType num_trials = 0;
290  IndexType num_threads = 1;
291  RandomNumberGenerator random_number_generator(num_threads);
292 
293  while (i < num_vectors)
294  {
295  if ((success == 0) && (num_trials >= max_num_trials))
296  {
297  std::cerr << "ERROR: Cannot orthogonalize vectors after " \
298  << num_trials << " trials. Aborting." \
299  << std::endl;
300  abort();
301  }
302 
303  // Reset on new trial (if it was set to 0 before to start a new trial)
304  success = 1;
305 
306  // j iterates on previous vectors in a window of at most vector_size
307  if (static_cast<LongIndexType>(i) > vector_size)
308  {
309  // When vector_size is smaller than i, it is fine to cast to signed
310  start_j = i - static_cast<IndexType>(vector_size);
311  }
312  else
313  {
314  start_j = 0;
315  }
316 
317  // Reorthogonalize against previous vectors
318  for (j=start_j; j < i; ++j)
319  {
320  // Norm of the j-th vector
322  &vectors[j*vector_size], vector_size);
323 
324  // Check norm
325  if (norm_j < epsilon * sqrt(vector_size))
326  {
327  std::cerr << "WARNING: norm of the given vector is too " \
328  << " small. Cannot reorthogonalize against zero" \
329  << "vector. Skipping."
330  << std::endl;
331  continue;
332  }
333 
334  // Projecting i-th vector to j-th vector
336  &vectors[i*vector_size], &vectors[j*vector_size],
337  vector_size);
338 
339  // Scale of subtraction
340  DataType scale = inner_prod / (norm_j * norm_j);
341 
342  // Subtraction
344  &vectors[vector_size*j], vector_size, scale,
345  &vectors[vector_size*i]);
346 
347  // Norm of the i-th vector
349  &vectors[i*vector_size], vector_size);
350 
351  // If the norm is too small, regenerate the i-th vector randomly
352  if (norm_i < epsilon * sqrt(vector_size))
353  {
354  // Regenerate new random vector for i-th vector
356  random_number_generator, &vectors[i*vector_size],
357  vector_size, num_threads);
358 
359  // Repeat the reorthogonalization for i-th vector against
360  // all previous vectors again.
361  success = 0;
362  ++num_trials;
363  break;
364  }
365  }
366 
367  if (success == 1)
368  {
369  ++i;
370 
371  // Reset if num_trials was incremented before.
372  num_trials = 0;
373  }
374  }
375 }
static void generate_random_array(RandomNumberGenerator &random_number_generator, DataType *array, const LongIndexType array_size, const IndexType num_threads)
Generates a pseudo-random array with Rademacher distribution where elements are either +1 or -1.
Generates 64-bit integers on multiple parallel threads.
int LongIndexType
Definition: types.h:60

References cVectorOperations< DataType >::euclidean_norm(), RandomArrayGenerator< DataType >::generate_random_array(), cVectorOperations< DataType >::inner_product(), and cVectorOperations< DataType >::subtract_scaled_vector().

Here is the call graph for this function:

The documentation for this class was generated from the following files: