zhetrd.c 15.5 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
/**
 *
 * @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 zhetrd.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
 * @author Gregoire Pichon
 * @author Mathieu Faverge
 * @date 2010-11-15
 * @precisions normal z -> s d c
 *
 **/
#include "control/common.h"
29 30 31
#if !defined(CHAMELEON_SIMULATION)
#include <coreblas/lapacke.h>
#endif
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
 *
 * @ingroup MORSE_Complex64_t
 *
 *  MORSE_zhetrd - reduces a complex Hermitian matrix A to real symmetric
 *  tridiagonal form S using a two-stage approach
 *  First stage: reduction to band tridiagonal form (unitary Q1);
 *  Second stage: reduction from band to tridiagonal form (unitary
 *  Q2).  Let Q = Q1 * Q2 be the global unitary transformation; Q**H *
 *  A * Q = S.
 *
 *******************************************************************************
 *
 * @param[in] jobz
 *          Intended usage:
 *          = MorseNoVec: computes tridiagonal only;
 *          = MorseVec: computes tridiagonal and generate the orthogonal matrix Q.
 *
 * @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] D
 *          On exit, the diagonal elements of the tridiagonal matrix:
 *          D(i) = A(i,i).
 *
 * @param[out] E
 *          On exit, he off-diagonal elements of the tridiagonal matrix:
 *          E(i) = A(i,i+1) if uplo = MorseUpper, E(i) = A(i+1,i) if uplo = MorseLower.
 *
 * @param[out] descT
 *          On entry, descriptor as return by MORSE_Alloc_Workspace_zhetrd
 *          On exit, contains auxiliary factorization data.
 *
 * @param[out] Q
 *          On exit, if jobz = MorseVec, then if return value = 0, Q
 *          contains the N-by-N unitary matrix Q.
 *          If jobz = MorseNoVec, then it is not referenced.
 *
 * @param[in] LDQ
 *          The leading dimension of the array Q. LDQ >= N.
 *
 *******************************************************************************
 *
 * @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_zhetrd_Tile
 * @sa MORSE_zhetrd_Tile_Async
 * @sa MORSE_chetrd
 * @sa MORSE_dsytrd
 * @sa MORSE_ssytrd
 *
 ******************************************************************************/
int MORSE_zhetrd(MORSE_enum jobz, MORSE_enum uplo, int N,
                 MORSE_Complex64_t *A, int LDA,
                 double *D,
                 double *E,
                 MORSE_desc_t *descT,
                 MORSE_Complex64_t *Q, int LDQ)
{
    int NB;
    int status;
    MORSE_context_t *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t request = MORSE_REQUEST_INITIALIZER;
Mathieu Faverge's avatar
Mathieu Faverge committed
128
    MORSE_desc_t descAl, descAt;
129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148

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

    /* Check input arguments */
    if (jobz != MorseNoVec && jobz != MorseVec) {
        morse_error("MORSE_zhetrd", "illegal value of jobz");
        return -1;
    }
    if (uplo != MorseLower && uplo != MorseUpper) {
        morse_error("MORSE_zhetrd", "illegal value of uplo");
        return -1;
    }
    if (N < 0) {
        morse_error("MORSE_zhetrd", "illegal value of N");
        return -2;
    }
149
    if (LDA < chameleon_max(1, N)) {
150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168
        morse_error("MORSE_zhetrd", "illegal value of LDA");
        return -4;
    }

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

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

    morse_sequence_create(morse, &sequence);

169
    /* Submit the matrix conversion */
170 171
    morse_zlap2tile( morse, &descAl, &descAt, MorseUpperLower,
                     A, NB, NB, LDA, N, N, N, sequence, &request );
172 173 174 175

    /* Call the tile interface */
    MORSE_zhetrd_Tile_Async(jobz, uplo, &descA, D, E, descT, Q, LDQ, sequence, &request);

Mathieu Faverge's avatar
Cleanup  
Mathieu Faverge committed
176
    /* Submit the matrix conversion back */
Mathieu Faverge's avatar
Mathieu Faverge committed
177 178 179
    morse_ztile2lap( morse, &descAl, &descAt,
                     MorseUpperLower, sequence, &request );

180
    morse_sequence_wait(morse, sequence);
Mathieu Faverge's avatar
Mathieu Faverge committed
181

Mathieu Faverge's avatar
Mathieu Faverge committed
182 183
    /* Cleanup the temporary data */
    morse_ztile2lap_cleanup( morse, &descAl, &descAt );
184 185 186 187 188

    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
}
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
 *
 * @ingroup MORSE_Complex64_t_Tile
 *
 *  MORSE_zhetrd_Tile - reduces a complex Hermitian matrix A to real symmetric
 *  tridiagonal form S using a two-stage approach
 *  First stage: reduction to band tridiagonal form (unitary Q1);
 *  Second stage: reduction from band to tridiagonal form (unitary Q2).
 *  Let Q = Q1 * Q2 be the global unitary transformation;
 *  Q**H * A * Q = S.
 *  Tile equivalent of MORSE_zhetrd().
 *  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 tridiagonal only;
 *          = MorseVec: computes tridiagonal and generate the orthogonal matrix Q.
 *
 * @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] D
 *          On exit, the diagonal elements of the tridiagonal matrix:
 *          D(i) = A(i,i).
 *
 * @param[out] E
 *          On exit, he off-diagonal elements of the tridiagonal matrix:
 *          E(i) = A(i,i+1) if uplo = MorseUpper,
 *          E(i) = A(i+1,i) if uplo = MorseLower.
 *
 * @param[out] T
 *          On exit, auxiliary factorization data.
 *
 * @param[out] Q
 *          On exit, if jobz = MorseVec, then if return value = 0, Q
 *          contains the N-by-N unitary matrix Q.
 *          If jobz = MorseNoVec, then it is not referenced.
 *
 * @param[in] LDQ
 *          The leading dimension of the array Q. LDQ >= N.
 *
 *******************************************************************************
 *
 * @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_zhetrd
 * @sa MORSE_zhetrd_Tile_Async
 * @sa MORSE_chetrd_Tile
 * @sa MORSE_dsytrd_Tile
 * @sa MORSE_ssytrd_Tile
 * @sa MORSE_zhetrd_Tile
 *
 ******************************************************************************/
int MORSE_zhetrd_Tile(MORSE_enum jobz, MORSE_enum uplo,
                      MORSE_desc_t *A, double *D, double *E,
                      MORSE_desc_t *T, MORSE_Complex64_t *Q, int LDQ)
{
    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_zhetrd_Tile", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    morse_sequence_create(morse, &sequence);
    MORSE_zhetrd_Tile_Async(jobz, uplo, A, D, E, T, Q, LDQ, sequence, &request);
Mathieu Faverge's avatar
Mathieu Faverge committed
286

287
    morse_sequence_wait(morse, sequence);
Mathieu Faverge's avatar
Mathieu Faverge committed
288

289 290 291 292 293
    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
}

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
 *
 * @ingroup MORSE_Complex64_t_Tile_Async
 *
 *  MORSE_zhetrd_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] 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_zhetrd
 * @sa MORSE_zhetrd_Tile
 * @sa MORSE_chetrd_Tile_Async
 * @sa MORSE_dsytrd_Tile_Async
 * @sa MORSE_ssytrd_Tile_Async
 *
 ******************************************************************************/
int MORSE_zhetrd_Tile_Async(MORSE_enum jobz,
                            MORSE_enum uplo,
                            MORSE_desc_t *A,
                            double *W,
                            double *E,
                            MORSE_desc_t *T,
                            MORSE_Complex64_t *Q, int LDQ,
                            MORSE_sequence_t *sequence, MORSE_request_t *request)
{
    MORSE_context_t *morse;
Mathieu Faverge's avatar
Mathieu Faverge committed
336 337
    MORSE_desc_t descAl, descAt;
    MORSE_desc_t descABl, descABt;
338 339
    int N, NB, LDAB;
    int status;
340
    MORSE_desc_t D, *Dptr = NULL;
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

    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zhetrd_Tile_Async", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    if (sequence == NULL) {
        morse_fatal_error("MORSE_zhetrd_Tile_Async", "NULL sequence");
        return MORSE_ERR_UNALLOCATED;
    }
    if (request == NULL) {
        morse_fatal_error("MORSE_zhetrd_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_zhetrd_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_zhetrd_Tile_Async", "invalid descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    /* Check input arguments */
    if (jobz != MorseNoVec && jobz != MorseVec) {
        morse_error("MORSE_zhetrd_Tile_Async", "illegal value of jobz");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (uplo != MorseLower && uplo != MorseUpper) {
        morse_error("MORSE_zhetrd_Tile_Async", "illegal value of uplo");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (descA.m != descA.n) {
        morse_error("MORSE_zhetrd_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_zhetrd_Tile_Async", "only square tiles supported");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    N  = descA.m;
    NB = descA.mb;
394 395 396 397 398 399
#if defined(CHAMELEON_COPY_DIAG)
    {
        morse_zdesc_alloc_diag(D, A->mb, A->nb, chameleon_min(A->m, A->n), A->nb, 0, 0, chameleon_min(A->m, A->n), A->nb, A->p, A->q);
        Dptr = &D;
    }
#endif
400
    /* Reduction to band. On exit, T contains reflectors */
401
    morse_pzhetrd_he2hb( uplo, A, T, Dptr,
402 403 404 405 406 407 408 409 410 411 412 413 414 415 416
                         sequence, request );

    LDAB = NB+1;

    /* Allocate band structure */
    morse_zdesc_alloc_diag( descAB,
                            LDAB, NB, /* mb, nb */
                            LDAB, N,  /* lm, ln */
                            0, 0,     /* i, j */
                            LDAB, N,  /* m, n */
                            1, 1 );

    /* Copy data into band structure */
    morse_pztile2band( uplo, A, &descAB,
                       sequence, request );
Mathieu Faverge's avatar
Mathieu Faverge committed
417

418 419 420
    morse_sequence_wait(morse, sequence);

    /* Reduce band matrix to tridiagonal matrix */
421
#if !defined(CHAMELEON_SIMULATION)
422 423 424 425 426 427 428 429 430
    status = LAPACKE_zhbtrd( LAPACK_COL_MAJOR,
                             morse_lapack_const(jobz),
                             morse_lapack_const(uplo),
                             N, NB,
                             (MORSE_Complex64_t *) descAB.mat, LDAB,
                             W, E, Q, LDQ );
    if (status != 0) {
        morse_error("MORSE_zhetrd_Tile_Async", "LAPACKE_zhbtrd failed");
    }
431
#endif /* !defined(CHAMELEON_SIMULATION) */
432
    if (Dptr != NULL) {
Mathieu Faverge's avatar
Mathieu Faverge committed
433
        morse_desc_mat_free( Dptr );
434
    }
Mathieu Faverge's avatar
Mathieu Faverge committed
435
    morse_ztile2lap_cleanup( morse, &descABl, &descABt );
Mathieu Faverge's avatar
Mathieu Faverge committed
436
    (void)D;
437 438
    return MORSE_SUCCESS;
}