zhemm.c 15 KB
Newer Older
1 2
/**
 *
3 4
 * @copyright (c) 2009-2014 The University of Tennessee and The University
 *                          of Tennessee Research Foundation.
5 6
 *                          All rights reserved.
 * @copyright (c) 2012-2014 Inria. All rights reserved.
7
 * @copyright (c) 2012-2014 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria, Univ. Bordeaux. All rights reserved.
8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28
 *
 **/

/**
 *
 * @file zhemm.c
 *
 *  MORSE computational routines
 *  MORSE is a software package provided by Univ. of Tennessee,
 *  Univ. of California Berkeley and Univ. of Colorado Denver
 *
 * @version 2.5.0
 * @comment This file has been automatically generated
 *          from Plasma 2.5.0 for MORSE 1.0.0
 * @author Mathieu Faverge
 * @author Emmanuel Agullo
 * @author Cedric Castagnede
 * @date 2010-11-15
 * @precisions normal z -> c
 *
 **/
29
#include "control/common.h"
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
 *
 * @ingroup MORSE_Complex64_t
 *
 *  MORSE_zhemm - Performs one of the matrix-matrix operations
 *
 *     \f[ C = \alpha \times A \times B + \beta \times C \f]
 *
 *  or
 *
 *     \f[ C = \alpha \times B \times A + \beta \times C \f]
 *
 *  where alpha and beta are scalars, A is an hermitian matrix and  B and
 *  C are m by n matrices.
 *
 *******************************************************************************
 *
 * @param[in] side
 *          Specifies whether the hermitian matrix A appears on the
 *          left or right in the operation as follows:
 *          = MorseLeft:      \f[ C = \alpha \times A \times B + \beta \times C \f]
 *          = MorseRight:     \f[ C = \alpha \times B \times A + \beta \times C \f]
 *
 * @param[in] uplo
 *          Specifies whether the upper or lower triangular part of
 *          the hermitian matrix A is to be referenced as follows:
 *          = MorseLower:     Only the lower triangular part of the
 *                             hermitian matrix A is to be referenced.
 *          = MorseUpper:     Only the upper triangular part of the
 *                             hermitian matrix A is to be referenced.
 *
 * @param[in] M
 *          Specifies the number of rows of the matrix C. M >= 0.
 *
 * @param[in] N
 *          Specifies the number of columns of the matrix C. N >= 0.
 *
 * @param[in] alpha
 *          Specifies the scalar alpha.
 *
 * @param[in] A
 *          A is a LDA-by-ka matrix, where ka is M when side = MorseLeft,
 *          and is N otherwise. Only the uplo triangular part is referenced.
 *
 * @param[in] LDA
 *          The leading dimension of the array A. LDA >= max(1,ka).
 *
 * @param[in] B
 *          B is a LDB-by-N matrix, where the leading M-by-N part of
 *          the array B must contain the matrix B.
 *
 * @param[in] LDB
 *          The leading dimension of the array B. LDB >= max(1,M).
 *
 * @param[in] beta
 *          Specifies the scalar beta.
 *
 * @param[in,out] C
 *          C is a LDC-by-N matrix.
 *          On exit, the array is overwritten by the M by N updated matrix.
 *
 * @param[in] LDC
 *          The leading dimension of the array C. LDC >= max(1,M).
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *
 *******************************************************************************
 *
 * @sa MORSE_zhemm_Tile
 * @sa MORSE_chemm
 * @sa MORSE_dhemm
 * @sa MORSE_shemm
 *
 ******************************************************************************/
int MORSE_zhemm(MORSE_enum side, MORSE_enum uplo, int M, int N,
110 111 112
                MORSE_Complex64_t alpha, MORSE_Complex64_t *A, int LDA,
                MORSE_Complex64_t *B, int LDB,
                MORSE_Complex64_t beta,  MORSE_Complex64_t *C, int LDC)
113 114 115 116
{
    int NB;
    int Am;
    int status;
117 118 119
    MORSE_desc_t descAl, descAt;
    MORSE_desc_t descBl, descBt;
    MORSE_desc_t descCl, descCt;
120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147
    MORSE_context_t *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t request = MORSE_REQUEST_INITIALIZER;

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

    /* Check input arguments */
    if ( (side != MorseLeft) && (side != MorseRight) ){
        morse_error("MORSE_zhemm", "illegal value of side");
        return -1;
    }
    if ((uplo != MorseLower) && (uplo != MorseUpper)) {
        morse_error("MORSE_zhemm", "illegal value of uplo");
        return -2;
    }
    Am = ( side == MorseLeft ) ? M : N;
    if (M < 0) {
        morse_error("MORSE_zhemm", "illegal value of M");
        return -3;
    }
    if (N < 0) {
        morse_error("MORSE_zhemm", "illegal value of N");
        return -4;
    }
148
    if (LDA < chameleon_max(1, Am)) {
149 150 151
        morse_error("MORSE_zhemm", "illegal value of LDA");
        return -7;
    }
152
    if (LDB < chameleon_max(1, M)) {
153 154 155
        morse_error("MORSE_zhemm", "illegal value of LDB");
        return -9;
    }
156
    if (LDC < chameleon_max(1, M)) {
157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177
        morse_error("MORSE_zhemm", "illegal value of LDC");
        return -12;
    }

    /* Quick return */
    if (M == 0 || N == 0 ||
        ((alpha == (MORSE_Complex64_t)0.0) && beta == (MORSE_Complex64_t)1.0))
        return MORSE_SUCCESS;

    /* Tune NB depending on M, N & NRHS; Set NBNBSIZE */
    status = morse_tune(MORSE_FUNC_ZHEMM, M, N, 0);
    if (status != MORSE_SUCCESS) {
        morse_error("MORSE_zhemm", "morse_tune() failed");
        return status;
    }

    /* Set MT & NT & KT */
    NB = MORSE_NB;

    morse_sequence_create(morse, &sequence);

178
    /* Submit the matrix conversion */
179 180 181 182 183 184
    morse_zlap2tile( morse, &descAl, &descAt, MorseUpperLower,
                     A, NB, NB, LDA, Am, Am, Am, sequence, &request );
    morse_zlap2tile( morse, &descBl, &descBt, MorseUpperLower,
                     B, NB, NB, LDB, N, M,  N, sequence, &request );
    morse_zlap2tile( morse, &descCl, &descCt, MorseUpperLower,
                     C, NB, NB, LDC, N, M,  N, sequence, &request );
185 186

    /* Call the tile interface */
187
    MORSE_zhemm_Tile_Async(  side, uplo, alpha, &descAt, &descBt, beta, &descCt, sequence, &request );
188

Mathieu Faverge's avatar
Mathieu Faverge committed
189
    /* Submit the matrix conversion back */
Mathieu Faverge's avatar
Mathieu Faverge committed
190 191 192
    morse_ztile2lap( morse, &descCl, &descCt,
                     MorseUpperLower, sequence, &request );

193
    morse_sequence_wait(morse, sequence);
Mathieu Faverge's avatar
Mathieu Faverge committed
194

Mathieu Faverge's avatar
Mathieu Faverge committed
195 196 197 198
    /* Cleanup the temporary data */
    morse_ztile2lap_cleanup( morse, &descAl, &descAt );
    morse_ztile2lap_cleanup( morse, &descBl, &descBt );
    morse_ztile2lap_cleanup( morse, &descCl, &descCt );
199 200 201 202 203 204

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

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
 *
 * @ingroup MORSE_Complex64_t_Tile
 *
 *  MORSE_zhemm_Tile - Performs Hermitian matrix multiplication.
 *  Tile equivalent of MORSE_zhemm().
 *  Operates on matrices stored by tiles.
 *  All matrices are passed through descriptors.
 *  All dimensions are taken from the descriptors.
 *
 *******************************************************************************
 *
 * @param[in] side
 *          Specifies whether the hermitian matrix A appears on the
 *          left or right in the operation as follows:
 *          = MorseLeft:      \f[ C = \alpha \times A \times B + \beta \times C \f]
 *          = MorseRight:     \f[ C = \alpha \times B \times A + \beta \times C \f]
 *
 * @param[in] uplo
 *          Specifies whether the upper or lower triangular part of
 *          the hermitian matrix A is to be referenced as follows:
 *          = MorseLower:     Only the lower triangular part of the
 *                             hermitian matrix A is to be referenced.
 *          = MorseUpper:     Only the upper triangular part of the
 *                             hermitian matrix A is to be referenced.
 *
 * @param[in] alpha
 *          Specifies the scalar alpha.
 *
 * @param[in] A
 *          A is a LDA-by-ka matrix, where ka is M when side = MorseLeft,
 *          and is N otherwise. Only the uplo triangular part is referenced.
 *
 * @param[in] B
 *          B is a LDB-by-N matrix, where the leading M-by-N part of
 *          the array B must contain the matrix B.
 *
 * @param[in] beta
 *          Specifies the scalar beta.
 *
 * @param[in,out] C
 *          C is a LDC-by-N matrix.
 *          On exit, the array is overwritten by the M by N updated matrix.
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *
 *******************************************************************************
 *
 * @sa MORSE_zhemm
 * @sa MORSE_zhemm_Tile_Async
 * @sa MORSE_chemm_Tile
 * @sa MORSE_dhemm_Tile
 * @sa MORSE_shemm_Tile
 *
 ******************************************************************************/
int MORSE_zhemm_Tile(MORSE_enum side, MORSE_enum uplo,
265 266
                     MORSE_Complex64_t alpha, MORSE_desc_t *A, MORSE_desc_t *B,
                     MORSE_Complex64_t beta,  MORSE_desc_t *C)
267 268 269 270 271 272 273 274 275 276 277 278
{
    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_zhemm_Tile", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    morse_sequence_create(morse, &sequence);
279
    MORSE_zhemm_Tile_Async(side, uplo, alpha, A, B, beta, C, sequence, &request );
280 281 282
    RUNTIME_desc_flush( A, sequence );
    RUNTIME_desc_flush( B, sequence );
    RUNTIME_desc_flush( C, sequence );
Mathieu Faverge's avatar
Mathieu Faverge committed
283

Mathieu Faverge's avatar
Mathieu Faverge committed
284
    morse_sequence_wait(morse, sequence);
285

286 287 288 289 290
    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
}

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
 *
 * @ingroup MORSE_Complex64_t_Tile_Async
 *
 *  MORSE_zhemm_Tile_Async - Performs Hermitian matrix multiplication.
 *  Non-blocking equivalent of MORSE_zhemm_Tile().
 *  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_zhemm
 * @sa MORSE_zhemm_Tile
 * @sa MORSE_chemm_Tile_Async
 * @sa MORSE_dhemm_Tile_Async
 * @sa MORSE_shemm_Tile_Async
 *
 ******************************************************************************/
int MORSE_zhemm_Tile_Async(MORSE_enum side, MORSE_enum uplo,
320 321 322
                           MORSE_Complex64_t alpha, MORSE_desc_t *A, MORSE_desc_t *B,
                           MORSE_Complex64_t beta,  MORSE_desc_t *C,
                           MORSE_sequence_t *sequence, MORSE_request_t *request)
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 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418
{
    MORSE_context_t *morse;

    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zhemm_Tile_Async", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    if (sequence == NULL) {
        morse_fatal_error("MORSE_zhemm_Tile_Async", "NULL sequence");
        return MORSE_ERR_UNALLOCATED;
    }
    if (request == NULL) {
        morse_fatal_error("MORSE_zhemm_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_zhemm_Tile_Async", "invalid first descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (morse_desc_check(B) != MORSE_SUCCESS) {
        morse_error("MORSE_zhemm_Tile_Async", "invalid second descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if (morse_desc_check(C) != MORSE_SUCCESS) {
        morse_error("MORSE_zhemm_Tile_Async", "invalid third descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    /* Check input arguments */
    if ( (side != MorseLeft) && (side != MorseRight) ){
        morse_error("MORSE_zhemm_Tile_Async", "illegal value of side");
        return morse_request_fail(sequence, request, -1);
    }
    if ((uplo != MorseLower) && (uplo != MorseUpper)) {
        morse_error("MORSE_zhemm_Tile_Async", "illegal value of uplo");
        return morse_request_fail(sequence, request, -2);
    }

    /* Check matrices sizes */
    if ( (B->m != C->m) || (B->n != C->n) ) {
        morse_error("MORSE_zhemm_Tile_Async", "B and C must have the same size");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if ( (A->m != A->n) ||
         ( (side == MorseLeft)  && (A->m != B->m ) ) ||
         ( (side == MorseRight) && (A->m != B->n ) ) ) {
        morse_error("MORSE_zhemm_Tile_Async", "Matrix A must be square of size M or N regarding side.");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    /* Check tiles sizes */
    if ( (B->mb != C->mb) || (B->nb != C->nb) ) {
        morse_error("MORSE_zhemm_Tile_Async", "B and C must have the same tile sizes");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
    if ( (A->mb != A->nb) ||
         ( (side == MorseLeft)  && (A->mb != B->mb ) ) ||
         ( (side == MorseRight) && (A->mb != B->nb ) ) ) {
        morse_error("MORSE_zhemm_Tile_Async", "Matrix A must be square with square tiles wich fits the reagding tile size of B and C");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    /* Check submatrix starting point */
    /* if ( (B->i != C->i) || (B->j != C->j) ) { */
    /*     morse_error("MORSE_zhemm_Tile_Async", "B and C submatrices doesn't match"); */
    /*     return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE); */
    /* } */
    /* if ( (A->i != A->j) ||  */
    /*          ( (side == MorseLeft)  && (A->i != B->i ) ) ||  */
    /*          ( (side == MorseRight) && (A->i != B->j ) ) ) { */
    /*     morse_error("MORSE_zhemm_Tile_Async", "Submatrix A must start on diagnonal and match submatrices B and C."); */
    /*     return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE); */
    /* } */
    if( (A->i != 0) || (A->j != 0) ||
        (B->i != 0) || (B->j != 0) ||
        (C->i != 0) || (C->j != 0) ) {
        morse_error("MORSE_zhemm_Tile_Async", "Submatrices are not supported for now");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }

    /* Quick return */
    if (C->m == 0 || C->n == 0 ||
        ( (alpha == (MORSE_Complex64_t)0.0) && (beta == (MORSE_Complex64_t)1.0) ))
        return MORSE_SUCCESS;

    morse_pzhemm(side, uplo, alpha, A, B, beta, C, sequence, request);

    return MORSE_SUCCESS;
}