imate
C++/CUDA Reference
Loading...
Searching...
No Matches
diagonalization.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// Headers
14// =======
15
16#include "./diagonalization.h"
17#include <cassert> // assert
18#include <cstdlib> // NULL
19#include "./lapack_api.h" // lapack_xstev, lapack_xbdsdc,
20
21
22// ================
23// eigh tridiagonal
24// ================
25
87
88template <typename DataType>
90 DataType* diagonals,
91 DataType* subdiagonals,
92 DataType* eigenvectors,
93 IndexType matrix_size)
94{
95 // 'V' computes both eigenvalues and eigenvectors
96 char jobz = 'V';
97
98 // Workspace array
99 DataType* work = new DataType[2*matrix_size - 2];
100
101 // Leading dimension of 2D eigenvectors array
102 int ldz = matrix_size;
103
104 // matrix size
105 int n = static_cast<int>(matrix_size);
106
107 // Error code output
108 int info;
109
110 // Calling Fortran subroutine
111 lapack_xstev(&jobz, &n, diagonals, subdiagonals, eigenvectors, &ldz, work,
112 &info);
113
114 delete[] work;
115 assert((info == 0, "?stev subroutine returned non-zero status."));
116
117 return info;
118}
119
120
121// ==============
122// svd bidiagonal
123// ==============
124
188
189template <typename DataType>
191 DataType* diagonals,
192 DataType* supdiagonals,
193 DataType* U,
194 DataType* Vt,
195 IndexType matrix_size)
196{
197 // Code 'U' indicates the matrix is upper bi-diagonal
198 char UPLO = 'U';
199
200 // Code 'I' computes both singular values and singular vectors
201 char COMPQ = 'I';
202
203 // matrix size
204 int n = static_cast<int>(matrix_size);
205
206 // Leading dimensions of arrays U and Vt
207 int LDU = matrix_size;
208 int LDVT = matrix_size;
209
210 // There arrays are not referenced when COMPQ is set to 'I'
211 DataType* Q = NULL;
212 int* IQ = NULL;
213
214 // Work arrays
215 DataType* work = new DataType[3*matrix_size*matrix_size + 4*matrix_size];
216 int* iwork = new int[8 * matrix_size];
217
218 // Error code output
219 int info;
220
221 // Calling Fortran subroutine
222 lapack_xbdsdc(&UPLO, &COMPQ, &n, diagonals, supdiagonals, U, &LDU, Vt,
223 &LDVT, Q, IQ, work, iwork, &info);
224
225 delete[] work;
226 delete[] iwork;
227
228 assert((info == 0, "?bdsdc subroutine returned non-zero status."));
229
230 return info;
231}
232
233
234// ===============================
235// Explicit template instantiation
236// ===============================
237
238template class Diagonalization<float>;
239template class Diagonalization<double>;
240template class Diagonalization<long double>;
A static class to find eigenvalues and eigenvectors (diagonalize) tridiagonal and bidiagonal matrices...
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...
void lapack_xbdsdc(char *uplo, char *compq, int *n, DataType *d, DataType *e, DataType *u, int *ldu, DataType *vt, int *ldvt, DataType *q, int *iq, DataType *work, int *iwork, int *info)
void lapack_xstev(char *jobz, int *n, DataType *d, DataType *e, DataType *z, int *ldz, DataType *work, int *info)
int IndexType
Definition types.h:65