imate
C++/CUDA Reference
c_orthogonalization.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 // Imports
14 // =======
15 
16 #include "./c_orthogonalization.h"
17 #include <cstdlib> // abort
18 #include <iostream> // std::cerr, std::endl
19 #include <cmath> // sqrt, std::fabs
20 #include <limits> // std::numeric_limits
21 #include "../_c_basic_algebra/c_vector_operations.h" // cVectorOperations
22 #include "../_random_generator/random_array_generator.h" // RandomArrayGene...
23 #include "../_random_generator/random_number_generator.h" // RandomNumberGe...
24 
25 
26 // ====================
27 // gram schmidt process
28 // ====================
29 
123 
124 template <typename DataType>
126  const DataType* V,
127  const LongIndexType vector_size,
128  const IndexType num_vectors,
129  const IndexType last_vector,
130  const FlagType num_ortho,
131  DataType* v)
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 }
225 
226 
227 // =====================
228 // orthogonalize vectors
229 // =====================
230 
267 
268 template <typename DataType>
270  DataType* vectors,
271  const LongIndexType vector_size,
272  const IndexType num_vectors)
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 }
376 
377 
378 // ===============================
379 // Explicit template instantiation
380 // ===============================
381 
382 template class cOrthogonalization<float>;
383 template class cOrthogonalization<double>;
384 template class cOrthogonalization<long double>;
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.
A static class for orthogonalization of vector bases. This class acts as a templated namespace,...
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 orthogonalize_vectors(DataType *vectors, const LongIndexType vector_size, const IndexType num_vectors)
Orthogonalizes set of vectors mutually using modified Gram-Schmidt process.
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 LongIndexType
Definition: types.h:60
int FlagType
Definition: types.h:68
int IndexType
Definition: types.h:65