zheevd.c 18.4 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
/**
 *
 * @copyright (c) 2009-2014 The University of Tennessee and The University
 *                          of Tennessee Research Foundation.
 *                          All rights reserved.
 * @copyright (c) 2012-2014 Inria. All rights reserved.
 * @copyright (c) 2012-2014 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria, Univ. Bordeaux. All rights reserved.
 *
 **/

/**
 * @file zheevd.c
 *
 *  PLASMA computational routines
 *  Release Date: November, 15th 2009
 *  PLASMA is a software package provided by Univ. of Tennessee,
 *  Univ. of California Berkeley and Univ. of Colorado Denver
 *
 * @version 2.7.1
 * @author Azzam Haidar
 * @author Hatem Ltaief
 * @date 2010-11-15
 * @precisions normal z -> s d c
 *
 **/
#include <coreblas/include/lapacke.h>
#include "control/common.h"

/***************************************************************************//**
 *
 * @ingroup MORSE_Complex64_t
 *
 *  MORSE_zheevd - Computes all eigenvalues and, optionally,
 *  eigenvectors of a complex Hermitian matrix A. The matrix A is
 *  preliminary reduced to tridiagonal form using a two-stage
 *  approach:
 *  First stage: reduction to band tridiagonal form;
 *  Second stage: reduction from band to tridiagonal form.
 *
 *******************************************************************************
 *
 * @param[in] jobz
 *          Intended usage:
 *          = MorseNoVec: computes eigenvalues only;
 *          = MorseVec: computes eigenvalues and eigenvectors.
 *
 * @param[in] uplo
 *          Specifies whether the matrix A is upper triangular or
 *          lower triangular:
 *          = MorseUpper: Upper triangle of A is stored;
 *          = MorseLower: Lower triangle of A is stored.
 *
 * @param[in] N
 *          The order of the matrix A. N >= 0.
 *
 * @param[in,out] A
 *          On entry, the symmetric (or Hermitian) matrix A.
 *          If uplo = MorseUpper, the leading N-by-N upper triangular
 *          part of A contains the upper triangular part of the matrix
 *          A, and the strictly lower triangular part of A is not
 *          referenced.
 *          If uplo = MorseLower, the leading N-by-N lower triangular
 *          part of A contains the lower triangular part of the matrix
 *          A, and the strictly upper triangular part of A is not
 *          referenced.
 *          On exit, the lower triangle (if uplo = MorseLower) or the
 *          upper triangle (if uplo = MorseUpper) of A, including the
 *          diagonal, is destroyed.
 *
 * @param[in] LDA
 *          The leading dimension of the array A. LDA >= max(1,N).
 *
 * @param[out] W
 *          On exit, if info = 0, the eigenvalues.
 *
 * @param[in, out] descT
 *          On entry, descriptor as return by MORSE_Alloc_Workspace_zheevd
 *          On exit, contains auxiliary factorization data.
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *          \retval <0 if -i, the i-th argument had an illegal value
 *          \retval >0 if INFO = i, the algorithm failed to converge; i
 *               off-diagonal elements of an intermediate tridiagonal
 *               form did not converge to zero.
 *
 *******************************************************************************
 *
 * @sa MORSE_zheevd_Tile
 * @sa MORSE_zheevd_Tile_Async
 * @sa MORSE_cheevd
 * @sa MORSE_dsyev
 * @sa MORSE_ssyev
 *
 ******************************************************************************/
int MORSE_zheevd(MORSE_enum jobz, MORSE_enum uplo, int N,
                 MORSE_Complex64_t *A, int LDA,
                 double *W,
                 MORSE_desc_t *descT)
{
    int NB;
    int status;
    MORSE_context_t  *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t   request = MORSE_REQUEST_INITIALIZER;
    MORSE_desc_t      descA;

    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zheevd", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }

    /* Check input arguments */
    if (jobz != MorseNoVec && jobz != MorseVec) {
        morse_error("MORSE_zheevd", "illegal value of jobz");
        return -1;
    }
    if (uplo != MorseLower && uplo != MorseUpper) {
        morse_error("MORSE_zheevd", "illegal value of uplo");
        return -2;
    }
    if (N < 0) {
        morse_error("MORSE_zheevd", "illegal value of N");
        return -3;
    }
129
    if (LDA < chameleon_max(1, N)) {
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
        morse_error("MORSE_zheevd", "illegal value of LDA");
        return -5;
    }

    /* Quick return */
    if (N == 0)
        return MORSE_SUCCESS;

    /* Tune NB & IB depending on N; Set NBNB */
    status = morse_tune(MORSE_FUNC_ZHEEVD, N, N, 0);
    if (status != MORSE_SUCCESS) {
        morse_error("MORSE_zheevd", "morse_tune() failed");
        return status;
    }

    /* Set NT */
    NB = MORSE_NB;

    morse_sequence_create(morse, &sequence);

    /* if ( MORSE_TRANSLATION == MORSE_OUTOFPLACE ) { */
    morse_zooplap2tile( descA, A, NB, NB, LDA, N, 0, 0, N, N, sequence, &request,
                        morse_desc_mat_free(&(descA)) );
    /* } else { */
    /*     morse_ziplap2tile( descA, A, NB, NB, LDA, N, 0, 0, N, N, */
    /*                         sequence, &requoest); */
    /* } */

    /* Call the tile interface */
    MORSE_zheevd_Tile_Async(jobz, uplo, &descA, W, descT, sequence, &request);

    /* if ( MORSE_TRANSLATION == MORSE_OUTOFPLACE ) { */
    morse_zooptile2lap( descA, A, NB, NB, LDA, N,  sequence, &request);
    morse_sequence_wait(morse, sequence);
    morse_desc_mat_free(&descA);
    /* } else { */
    /*     morse_ziptile2lap( descA, A, NB, NB, LDA, N,  sequence, &request); */
    /*     morse_sequence_wait(morse, sequence); */
    /* } */

    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
}
/***************************************************************************//**
 *
 * @ingroup MORSE_Complex64_t_Tile
 *
 *  MORSE_zheevd_Tile - Computes all eigenvalues and, optionally, eigenvectors of a
 *  complex Hermitian matrix A using a two-stage approach:
 *  First stage: reduction to band tridiagonal form;
 *  Second stage: reduction from band to tridiagonal form.
 *
 *  Operates on matrices stored by tiles.
 *  All matrices are passed through descriptors.
 *  All dimensions are taken from the descriptors.
 *
 *******************************************************************************
 *
 * @param[in] jobz
 *          Intended usage:
 *          = MorseNoVec: computes eigenvalues only;
 *          = MorseVec: computes eigenvalues and eigenvectors.
 *
 * @param[in] uplo
 *          Specifies whether the matrix A is upper triangular or
 *          lower triangular:
 *          = MorseUpper: Upper triangle of A is stored;
 *          = MorseLower: Lower triangle of A is stored.
 *
 * @param[in,out] A
 *          On entry, the symmetric (or Hermitian) matrix A.
 *          If uplo = MorseUpper, the leading N-by-N upper triangular
 *          part of A contains the upper triangular part of the matrix
 *          A, and the strictly lower triangular part of A is not
 *          referenced.
 *          If UPLO = 'L', the leading N-by-N lower triangular part of
 *          A contains the lower triangular part of the matrix A, and
 *          the strictly upper triangular part of A is not referenced.
 *          On exit, if jobz = MorseVec, then if return value = 0, A
 *          contains the orthonormal eigenvectors of the matrix A.
 *          If jobz = MorseNoVec, then on exit the lower triangle (if
 *          uplo = MorseLower) or the upper triangle (if uplo =
 *          MorseUpper) of A, including the diagonal, is destroyed.*
 *
 * @param[out] W
 *          On exit, if info = 0, the eigenvalues.
 *
 * @param[in,out] T
 *          On entry, descriptor as return by
 *          MORSE_Alloc_Workspace_zheevd
 *          On exit, contains auxiliary factorization data.
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *          \retval <0 if -i, the i-th argument had an illegal value
 *          \retval >0 if INFO = i, the algorithm failed to converge; i
 *               off-diagonal elements of an intermediate tridiagonal
 *               form did not converge to zero.
 *
 *******************************************************************************
 *
 * @sa MORSE_zheevd_Tile
 * @sa MORSE_zheevd_Tile_Async
 * @sa MORSE_cheevd
 * @sa MORSE_dsyev
 * @sa MORSE_ssyev
 *
 ******************************************************************************/
int MORSE_zheevd_Tile(MORSE_enum jobz, MORSE_enum uplo,
                      MORSE_desc_t *A, double *W, MORSE_desc_t *T)
{
    MORSE_context_t  *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t   request = MORSE_REQUEST_INITIALIZER;
    int status;

    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zheevd_Tile", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    morse_sequence_create(morse, &sequence);
    MORSE_zheevd_Tile_Async(jobz, uplo, A, W, T, sequence, &request);
    morse_sequence_wait(morse, sequence);

    RUNTIME_desc_getoncpu(A);
    RUNTIME_desc_getoncpu(T);

    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
}

/***************************************************************************//**
 *
 * @ingroup MORSE_Complex64_t_Tile_Async
 *
 *  MORSE_zheevd_Tile_Async - Computes all eigenvalues and,
 *  optionally, eigenvectors of a complex Hermitian matrix A using a
 *  two-stage approach:
 *  First stage: reduction to band tridiagonal form;
 *  Second stage: reduction from band to tridiagonal form.
 *
 *  May return before the computation is finished.
 *  Allows for pipelining of operations at runtime.
 *
 *******************************************************************************
 *
 * @param[in] jobz
 *          Intended usage:
 *          = MorseNoVec: computes eigenvalues only;
 *          = MorseVec: computes eigenvalues and eigenvectors.
 *
 * @param[in] uplo
 *          Specifies whether the matrix A is upper triangular or
 *          lower triangular:
 *          = MorseUpper: Upper triangle of A is stored;
 *          = MorseLower: Lower triangle of A is stored.
 *
 * @param[in,out] A
 *          On entry, the symmetric (or Hermitian) matrix A.
 *          If uplo = MorseUpper, the leading N-by-N upper triangular
 *          part of A contains the upper triangular part of the matrix
 *          A, and the strictly lower triangular part of A is not
 *          referenced.
 *          If UPLO = 'L', the leading N-by-N lower triangular part of
 *          A contains the lower triangular part of the matrix A, and
 *          the strictly upper triangular part of A is not referenced.
 *          On exit, if jobz = MorseVec, then if return value = 0, A
 *          contains the orthonormal eigenvectors of the matrix A.
 *          If jobz = MorseNoVec, then on exit the lower triangle (if
 *          uplo = MorseLower) or the upper triangle (if uplo =
 *          MorseUpper) of A, including the diagonal, is destroyed.*
 *
 * @param[out] W
 *          On exit, if info = 0, the eigenvalues.
 *
 * @param[in,out] T
 *          On entry, descriptor as return by
 *          MORSE_Alloc_Workspace_zheevd
 *          On exit, contains auxiliary factorization data.
 *
 * @param[in] sequence
 *          Identifies the sequence of function calls that this call belongs to
 *          (for completion checks and exception handling purposes).
 *
 * @param[out] request
 *          Identifies this function call (for exception handling purposes).
 *
 *******************************************************************************
 *
 * @sa MORSE_zheevd
 * @sa MORSE_zheevd_Tile
 * @sa MORSE_cheevd_Tile_Async
 * @sa MORSE_dsyev_Tile_Async
 * @sa MORSE_ssyev_Tile_Async
 *
 ******************************************************************************/
int MORSE_zheevd_Tile_Async(MORSE_enum jobz, MORSE_enum uplo,
                            MORSE_desc_t *A, double *W, MORSE_desc_t *T,
                            MORSE_sequence_t *sequence, MORSE_request_t *request)
{
    MORSE_context_t *morse;
    MORSE_desc_t descA;
    MORSE_desc_t descT;
    MORSE_Complex64_t *Q2;
    int N, NB, status;
    double *E;
    MORSE_Complex64_t *V;
    MORSE_desc_t descQ2, descV;
    MORSE_desc_t *subA, *subQ, *subT;

    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zheevd_Tile_Async", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    if (sequence == NULL) {
        morse_fatal_error("MORSE_zheevd_Tile_Async", "NULL sequence");
        return MORSE_ERR_UNALLOCATED;
    }
    if (request == NULL) {
        morse_fatal_error("MORSE_zheevd_Tile_Async", "NULL request");
        return MORSE_ERR_UNALLOCATED;
    }
    /* Check sequence status */
    if (sequence->status == MORSE_SUCCESS)
        request->status = MORSE_SUCCESS;
    else
        return morse_request_fail(sequence, request, MORSE_ERR_SEQUENCE_FLUSHED);

    /* Check descriptors for correctness */
    if (morse_desc_check(A) != MORSE_SUCCESS) {
        morse_error("MORSE_zheevd_Tile_Async", "invalid descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    } else {
        descA = *A;
    }
    if (morse_desc_check(T) != MORSE_SUCCESS) {
        morse_error("MORSE_zheevd_Tile_Async", "invalid descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    } else {
        descT = *T;
    }
    /* Check input arguments */
    if (jobz != MorseNoVec && jobz != MorseVec) {
        morse_error("MORSE_zheevd_Tile_Async", "illegal value of jobz");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (uplo != MorseLower && uplo != MorseUpper) {
        morse_error("MORSE_zheevd_Tile_Async", "illegal value of uplo");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (descA.m != descA.n) {
        morse_error("MORSE_zheevd_Tile_Async", "matrix need to be square");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (descA.nb != descA.mb) {
        morse_error("MORSE_zheevd_Tile_Async", "only square tiles supported");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    N   = descA.m;
396
    NB  = chameleon_min(descA.mb,descA.m);
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509

    /* Allocate data structures for reduction to tridiagonal form */
    E = malloc( (N - 1) * sizeof(double) );
    if (E == NULL) {
        morse_error("MORSE_zheevd_Tile_Async", "malloc(E) failed");
        free(E);
        return MORSE_ERR_OUT_OF_RESOURCES;
    }

    if (jobz == MorseVec){
        /* Have to synchrone right now */
        Q2 = malloc( N * N * sizeof(MORSE_Complex64_t));
        /* For bug in lapacke */
        memset( Q2, 0, N * N * sizeof(MORSE_Complex64_t));
    }

    status = MORSE_zhetrd_Tile_Async( jobz, uplo,
                                      A, W, E, T,
                                      Q2, N,
                                      sequence, request );
    if (status != 0) {
        morse_error("MORSE_zheevd_Tile", "MORSE_zhetrd failed");
    }


    if (jobz == MorseNoVec){
        /* Tridiagonal eigensolver */
        status = LAPACKE_zstedc( LAPACK_COL_MAJOR,
                                 morse_lapack_const(jobz),
                                 N, W, E,
                                 NULL, N );
        if (status != 0) {
            morse_error("MORSE_zheevd_Tile_Async", "LAPACKE_zstedc failed");
        }
        free(E);
        return MORSE_SUCCESS;
    }

    V = malloc( N * N * sizeof(MORSE_Complex64_t) );
    if (V == NULL) {
        morse_error("MORSE_zheevd_Tile_Async", "malloc(V) failed");
        free(V);
        return MORSE_ERR_OUT_OF_RESOURCES;
    }
    /* For bug in lapacke */
    memset(V, 0, N * N * sizeof(MORSE_Complex64_t));

    /* Tridiagonal eigensolver */
    /* V contains the eigenvectors of the tridiagonal system on exit */
    status = LAPACKE_zstedc( LAPACK_COL_MAJOR,
                             morse_lapack_const(MorseIvec),
                             N, W, E,
                             V, N );
    if (status != 0) {
        morse_error("MORSE_zheevd_Tile_Async", "LAPACKE_zstedc failed");
    }

    /* Convert matrices in tile format */
    /* A/T from MORSE_zhetrd   refer  to Q1 (tile   layout) */
    /* Q   from MORSE_zhetrd   refers to Q2 (lapack layout) */
    /* V   from LAPACKE_zstedc refers to V  (lapack layout) */
    /* The final eigenvectors are (Q1 Q2 V) or (Q1^h Q2 V)  */
    morse_zooplap2tile( descQ2, Q2, NB, NB, N, N, 0, 0, N, N, sequence, request,
                        morse_desc_mat_free(&(descQ2)) );
    morse_zooplap2tile( descV,  V,  NB, NB, N, N, 0, 0, N, N, sequence, request,
                        morse_desc_mat_free(&(descQ2)); morse_desc_mat_free(&(descV)) );
    if (uplo == MorseLower)
    {
        subA = morse_desc_submatrix(&descA,  descA.mb,  0, descA.m -descA.mb,  descA.n-descA.nb);
        subQ = morse_desc_submatrix(&descQ2, descQ2.mb, 0, descQ2.m-descQ2.mb, descQ2.n        );
        subT = morse_desc_submatrix(&descT,  descT.mb,  0, descT.m -descT.mb,  descT.n-descT.nb);

        /* Compute Q2 = Q1 * Q2 */
        morse_pzunmqr( MorseLeft, MorseNoTrans,
                       subA, subQ, subT,
                       sequence, request );

        /* Compute the final eigenvectors A = (Q1 * Q2) * V */
        morse_pzgemm( MorseNoTrans, MorseNoTrans,
                      1.0, &descQ2, &descV,
                      0.0, &descA,
                      sequence, request);

    }
    else {
        subA = morse_desc_submatrix(&descA,  0, descA.nb,  descA.m -descA.mb,  descA.n -descA.nb );
        subQ = morse_desc_submatrix(&descQ2, descQ2.mb, 0, descQ2.m-descQ2.mb, descQ2.n          );
        subT = morse_desc_submatrix(&descT,  0, descT.nb,  descT.m -descT.mb,  descT.n -descT.nb );

        /* Compute Q2 = Q1^h * Q2 */
        morse_pzunmlq( MorseLeft, MorseConjTrans,
                       subA, subQ, subT,
                       sequence, request );

        /* Compute the final eigenvectors A =  (Q1^h * Q2) * V */
        morse_pzgemm( MorseNoTrans, MorseNoTrans,
                      1.0, &descQ2, &descV,
                      0.0, &descA,
                      sequence, request);
    }

    morse_sequence_wait(morse, sequence);

    free(subA); free(subQ); free(subT);
    morse_desc_mat_free(&descQ2);
    free(Q2);

    morse_desc_mat_free(&descV);
    free(V);

    free(E);
    return MORSE_SUCCESS;
}