21 #include "../_cu_basic_algebra/cu_vector_operations.h"
22 #include "../_random_generator/random_array_generator.h"
23 #include "../_random_generator/random_number_generator.h"
24 #include "../_cuda_utilities/cuda_interface.h"
127 template <
typename DataType>
129 cublasHandle_t cublas_handle,
139 if ((num_ortho == 0) || (num_vectors < 2))
144 else if ((num_ortho < 0) ||
145 (num_ortho >
static_cast<FlagType>(num_vectors)))
148 num_steps = num_vectors;
153 num_steps = num_ortho;
159 if (num_steps >
static_cast<IndexType>(vector_size))
161 num_steps = vector_size;
168 DataType epsilon = std::numeric_limits<DataType>::epsilon();
172 for (
IndexType step=0; step < num_steps; ++step)
175 if ((last_vector % num_vectors) >= step)
177 i = (last_vector % num_vectors) - step;
182 i = (last_vector % num_vectors) - step + num_vectors;
187 cublas_handle, &V[vector_size*i], vector_size);
190 if (norm < epsilon * sqrt(vector_size))
192 std::cerr <<
"WARNING: norm of the given vector is too small. " \
193 <<
"Cannot orthogonalize against zero vector. " \
194 <<
"Skipping." << std::endl;
200 cublas_handle, &V[vector_size*i], v, vector_size);
203 DataType scale = inner_prod / (norm * norm);
207 if (std::abs(scale - 1.0) <= 2.0 * epsilon)
211 cublas_handle, v, vector_size);
214 distance2 = norm_v*norm_v - 2.0*inner_prod + norm*norm;
218 if (distance2 < 2.0 * epsilon * vector_size)
226 cublas_handle, &V[vector_size*i], vector_size, scale, v);
274 template <
typename DataType>
276 cublasHandle_t cublas_handle,
293 DataType epsilon = std::numeric_limits<DataType>::epsilon();
299 DataType* buffer = NULL;
301 while (i < num_vectors)
303 if ((success == 0) && (num_trials >= max_num_trials))
305 std::cerr <<
"ERROR: Cannot orthogonalize vectors after " \
306 << num_trials <<
" trials. Aborting." \
318 start_j = i -
static_cast<IndexType>(vector_size);
326 for (j=start_j; j < i; ++j)
330 cublas_handle, &vectors[j*vector_size], vector_size);
333 if (norm_j < epsilon * sqrt(vector_size))
335 std::cerr <<
"WARNING: norm of the given vector is too " \
336 <<
" small. Cannot reorthogonalize against zero" \
337 <<
"vector. Skipping."
344 cublas_handle, &vectors[i*vector_size],
345 &vectors[j*vector_size], vector_size);
348 DataType scale = inner_prod / (norm_j * norm_j);
352 cublas_handle, &vectors[vector_size*j], vector_size, scale,
353 &vectors[vector_size*i]);
357 cublas_handle, &vectors[i*vector_size], vector_size);
360 if (norm_i < epsilon * sqrt(vector_size))
365 buffer =
new DataType[vector_size];
370 random_number_generator, buffer,
371 vector_size, num_threads);
375 buffer, vector_size, &vectors[i*vector_size]);
static void copy_to_device(const ArrayType *host_array, const LongIndexType array_size, ArrayType *device_array)
Copies memory on host to device memory.
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(cublasHandle_t cublas_handle, 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(cublasHandle_t cublas_handle, 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(cublasHandle_t cublas_handle, 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(cublasHandle_t cublas_handle, const DataType *vector1, const DataType *vector2, const LongIndexType vector_size)
Computes Euclidean inner product of two vectors.
static DataType euclidean_norm(cublasHandle_t cublas_handle, const DataType *vector, const LongIndexType vector_size)
Computes the Euclidean 2-norm of a 1D array.