22 #include "../_random_generator/random_array_generator.h"
25 #include "../_utilities/timer.h"
191 template <
typename DataType>
194 DataType* parameters,
198 const DataType exponent,
202 const DataType lanczos_tol,
205 const DataType error_atol,
206 const DataType error_rtol,
207 const DataType confidence_level,
208 const DataType outlier_significance_level,
217 float& alg_wall_time)
223 omp_set_num_threads(num_threads);
230 IndexType random_vectors_size = matrix_size * num_threads;
231 DataType* random_vectors =
new DataType[random_vectors_size];
246 unsigned int chunk_size =
static_cast<int>(
247 sqrt(
static_cast<DataType
>(max_num_samples) / num_threads));
259 #pragma omp parallel for schedule(dynamic, chunk_size)
260 for (i=0; i < max_num_samples; ++i)
262 if (!
static_cast<bool>(all_converged))
264 int thread_id = omp_get_thread_num();
268 A, parameters, num_inquiries, matrix_function, gram,
269 exponent, orthogonalize, lanczos_degree, lanczos_tol,
270 random_number_generator,
271 &random_vectors[matrix_size*thread_id], converged,
278 processed_samples_indices[num_processed_samples] = i;
279 ++num_processed_samples;
285 samples, min_num_samples, num_inquiries,
286 processed_samples_indices, num_processed_samples,
287 confidence_level, error_atol, error_rtol, error,
288 num_samples_used, converged);
295 alg_wall_time = timer.
elapsed();
299 confidence_level, outlier_significance_level, num_inquiries,
300 max_num_samples, num_samples_used, processed_samples_indices,
301 samples, num_outliers, trace, error);
304 delete[] random_vectors;
306 return all_converged;
410 template <
typename DataType>
413 DataType* parameters,
417 const DataType exponent,
420 const DataType lanczos_tol,
422 DataType* random_vector,
424 DataType* trace_estimate)
435 random_number_generator, random_vector, matrix_size, num_threads);
438 DataType* alpha =
new DataType[lanczos_degree];
439 DataType* beta =
new DataType[lanczos_degree];
443 DataType* eigenvectors = NULL;
444 DataType* left_singularvectors = NULL;
445 DataType* right_singularvectors_transposed = NULL;
448 IndexType required_num_inquiries = num_inquiries;
456 required_num_inquiries = 1;
462 DataType** theta =
new DataType*[num_inquiries];
463 for (j=0; j < num_inquiries; ++j)
465 theta[j] =
new DataType[lanczos_degree];
468 for (i=0; i < lanczos_degree; ++i)
475 DataType** tau =
new DataType*[num_inquiries];
476 for (j=0; j < num_inquiries; ++j)
478 tau[j] =
new DataType[lanczos_degree];
481 for (i=0; i < lanczos_degree; ++i)
497 for (j=0; j < required_num_inquiries; ++j)
503 if ((converged[j] == 1) && (required_num_inquiries == num_inquiries))
515 A, random_vector, matrix_size, lanczos_degree, lanczos_tol,
516 orthogonalize, alpha, beta);
519 left_singularvectors = \
520 new DataType[lanczos_size[j] * lanczos_size[j]];
521 right_singularvectors_transposed = \
522 new DataType[lanczos_size[j] * lanczos_size[j]];
526 alpha, beta, left_singularvectors,
527 right_singularvectors_transposed, lanczos_size[j]);
530 for (i=0; i < lanczos_size[j]; ++i)
532 theta[j][i] = alpha[i] * alpha[i];
533 tau[j][i] = right_singularvectors_transposed[i];
540 A, random_vector, matrix_size, lanczos_degree, lanczos_tol,
541 orthogonalize, alpha, beta);
544 eigenvectors =
new DataType[lanczos_size[j] * lanczos_size[j]];
548 alpha, beta, eigenvectors, lanczos_size[j]);
551 for (i=0; i < lanczos_size[j]; ++i)
553 theta[j][i] = alpha[i];
554 tau[j][i] = eigenvectors[i * lanczos_size[j]];
571 for (j=0; j < num_inquiries; ++j)
577 for (j=1; j < num_inquiries; ++j)
580 lanczos_size[j] = lanczos_size[0];
582 for (i=0; i < lanczos_size[j]; ++i)
588 ¶meters[j*num_parameters]);
591 tau[j][i] = tau[0][i];
597 DataType quadrature_sum;
598 for (j=0; j < num_inquiries; ++j)
601 if (converged[j] == 1)
607 quadrature_sum = 0.0;
614 for (i=0; i < lanczos_size[j]; ++i)
616 quadrature_sum += tau[j][i] * tau[j][i] * \
617 matrix_function->function(pow(theta[j][i], exponent));
620 trace_estimate[j] = matrix_size * quadrature_sum;
626 delete[] lanczos_size;
628 for (j=0; j < required_num_inquiries; ++j)
634 for (j=0; j < required_num_inquiries; ++j)
640 if (eigenvectors != NULL)
642 delete[] eigenvectors;
645 if (left_singularvectors != NULL)
647 delete[] left_singularvectors;
650 if (right_singularvectors_transposed != NULL)
652 delete[] right_singularvectors_transposed;
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.
IndexType c_lanczos_tridiagonalization(cLinearOperator< DataType > *A, const DataType *v, const LongIndexType n, const IndexType m, const DataType lanczos_tol, const FlagType orthogonalize, DataType *alpha, DataType *beta)
Tri-diagonalizes matrix A to T using the start vector v. is the Lanczos degree, which will be the siz...
static int eigh_tridiagonal(DataType *diagonals, DataType *subdiagonals, DataType *eigenvectors, IndexType matrix_size)
Computes all eigenvalues and eigenvectors of a real and symmetric tridiagonal matrix.
static int svd_bidiagonal(DataType *diagonals, DataType *subdiagonals, DataType *U, DataType *Vt, IndexType matrix_size)
Computes all singular-values and left and right eigenvectors of a real and symmetric upper bi-diagona...
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.
Records elasped wall time between two events.
void start()
Starts the timer.
void stop()
Stops the timer.
double elapsed() const
Returns the elapsed time in seconds.
Base class for linear operators. This class serves as interface for all derived classes.
virtual DataType get_eigenvalue(const DataType *known_parameters, const DataType known_eigenvalue, const DataType *inquiry_parameters) const =0
void set_parameters(DataType *parameters_)
Sets the scalar parameter this->parameters. Parameter is initialized to NULL. However,...
LongIndexType get_num_rows() const
FlagType is_eigenvalue_relation_known() const
Returns a flag that determines whether a relation between the parameters of the operator and its eige...
IndexType get_num_parameters() const
A static class to compute the trace of implicit matrix functions using stochastic Lanczos quadrature ...
static void _c_stochastic_lanczos_quadrature(cLinearOperator< DataType > *A, DataType *parameters, const IndexType num_inquiries, const Function *matrix_function, const FlagType gram, const DataType exponent, const FlagType orthogonalize, const IndexType lanczos_degree, const DataType lanczos_tol, RandomNumberGenerator &random_number_generator, DataType *random_vector, FlagType *converged, DataType *trace_estimate)
For a given random input vector, computes one Monte-Carlo sample to estimate trace using Lanczos quad...
static FlagType c_trace_estimator(cLinearOperator< DataType > *A, DataType *parameters, const IndexType num_inquiries, const Function *matrix_function, const FlagType gram, const DataType exponent, const FlagType orthogonalize, const int64_t seed, const IndexType lanczos_degree, const DataType lanczos_tol, const IndexType min_num_samples, const IndexType max_num_samples, const DataType error_atol, const DataType error_rtol, const DataType confidence_level, const DataType outlier_significance_level, const IndexType num_threads, DataType *trace, DataType *error, DataType **samples, IndexType *processed_samples_indices, IndexType *num_samples_used, IndexType *num_outliers, FlagType *converged, float &alg_wall_time)
Stochastic Lanczos quadrature method to estimate trace of a function of a linear operator....