zhemm.c 15.5 KB
Newer Older
1
/**
2 3
 *
 * @file zhemm.c
4
 *
Mathieu Faverge's avatar
Mathieu Faverge committed
5 6
 * @copyright 2009-2014 The University of Tennessee and The University of
 *                      Tennessee Research Foundation. All rights reserved.
Mathieu Faverge's avatar
Mathieu Faverge committed
7
 * @copyright 2012-2018 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8
 *                      Univ. Bordeaux. All rights reserved.
9
 *
10
 ***
11
 *
12
 * @brief Chameleon zhemm wrappers
13
 *
Mathieu Faverge's avatar
Mathieu Faverge committed
14
 * @version 1.0.0
15
 * @comment This file has been automatically generated
16
 *          from Plasma 2.5.0 for CHAMELEON 1.0.0
17 18 19 20 21 22
 * @author Mathieu Faverge
 * @author Emmanuel Agullo
 * @author Cedric Castagnede
 * @date 2010-11-15
 * @precisions normal z -> c
 *
23
 */
24
#include "control/common.h"
25

26 27
/**
 ********************************************************************************
28
 *
29
 * @ingroup CHAMELEON_Complex64_t
30
 *
31
 *  CHAMELEON_zhemm - Performs one of the matrix-matrix operations
32 33 34 35 36 37 38 39 40 41 42 43 44 45 46
 *
 *     \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:
47 48
 *          = ChamLeft:      \f[ C = \alpha \times A \times B + \beta \times C \f]
 *          = ChamRight:     \f[ C = \alpha \times B \times A + \beta \times C \f]
49 50 51 52
 *
 * @param[in] uplo
 *          Specifies whether the upper or lower triangular part of
 *          the hermitian matrix A is to be referenced as follows:
53
 *          = ChamLower:     Only the lower triangular part of the
54
 *                             hermitian matrix A is to be referenced.
55
 *          = ChamUpper:     Only the upper triangular part of the
56 57 58 59 60 61 62 63 64 65 66 67
 *                             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
68
 *          A is a LDA-by-ka matrix, where ka is M when side = ChamLeft,
69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92
 *          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).
 *
 *******************************************************************************
 *
93
 * @retval CHAMELEON_SUCCESS successful exit
94 95 96
 *
 *******************************************************************************
 *
97 98 99 100
 * @sa CHAMELEON_zhemm_Tile
 * @sa CHAMELEON_chemm
 * @sa CHAMELEON_dhemm
 * @sa CHAMELEON_shemm
101
 *
102
 */
103 104 105 106
int CHAMELEON_zhemm( cham_side_t side, cham_uplo_t uplo, int M, int N,
                 CHAMELEON_Complex64_t alpha, CHAMELEON_Complex64_t *A, int LDA,
                 CHAMELEON_Complex64_t *B, int LDB,
                 CHAMELEON_Complex64_t beta,  CHAMELEON_Complex64_t *C, int LDC )
107 108 109 110
{
    int NB;
    int Am;
    int status;
111 112 113
    CHAM_desc_t descAl, descAt;
    CHAM_desc_t descBl, descBt;
    CHAM_desc_t descCl, descCt;
Mathieu Faverge's avatar
Mathieu Faverge committed
114
    CHAM_context_t *chamctxt;
115 116
    RUNTIME_sequence_t *sequence = NULL;
    RUNTIME_request_t request = RUNTIME_REQUEST_INITIALIZER;
117

Mathieu Faverge's avatar
Mathieu Faverge committed
118 119 120
    chamctxt = chameleon_context_self();
    if (chamctxt == NULL) {
        chameleon_fatal_error("CHAMELEON_zhemm", "CHAMELEON not initialized");
121
        return CHAMELEON_ERR_NOT_INITIALIZED;
122 123 124
    }

    /* Check input arguments */
125
    if ( (side != ChamLeft) && (side != ChamRight) ){
Mathieu Faverge's avatar
Mathieu Faverge committed
126
        chameleon_error("CHAMELEON_zhemm", "illegal value of side");
127 128
        return -1;
    }
129
    if ((uplo != ChamLower) && (uplo != ChamUpper)) {
Mathieu Faverge's avatar
Mathieu Faverge committed
130
        chameleon_error("CHAMELEON_zhemm", "illegal value of uplo");
131 132
        return -2;
    }
133
    Am = ( side == ChamLeft ) ? M : N;
134
    if (M < 0) {
Mathieu Faverge's avatar
Mathieu Faverge committed
135
        chameleon_error("CHAMELEON_zhemm", "illegal value of M");
136 137 138
        return -3;
    }
    if (N < 0) {
Mathieu Faverge's avatar
Mathieu Faverge committed
139
        chameleon_error("CHAMELEON_zhemm", "illegal value of N");
140 141
        return -4;
    }
142
    if (LDA < chameleon_max(1, Am)) {
Mathieu Faverge's avatar
Mathieu Faverge committed
143
        chameleon_error("CHAMELEON_zhemm", "illegal value of LDA");
144 145
        return -7;
    }
146
    if (LDB < chameleon_max(1, M)) {
Mathieu Faverge's avatar
Mathieu Faverge committed
147
        chameleon_error("CHAMELEON_zhemm", "illegal value of LDB");
148 149
        return -9;
    }
150
    if (LDC < chameleon_max(1, M)) {
Mathieu Faverge's avatar
Mathieu Faverge committed
151
        chameleon_error("CHAMELEON_zhemm", "illegal value of LDC");
152 153 154 155 156
        return -12;
    }

    /* Quick return */
    if (M == 0 || N == 0 ||
157 158
        ((alpha == (CHAMELEON_Complex64_t)0.0) && beta == (CHAMELEON_Complex64_t)1.0))
        return CHAMELEON_SUCCESS;
159 160

    /* Tune NB depending on M, N & NRHS; Set NBNBSIZE */
Mathieu Faverge's avatar
Mathieu Faverge committed
161
    status = chameleon_tune(CHAMELEON_FUNC_ZHEMM, M, N, 0);
162
    if (status != CHAMELEON_SUCCESS) {
Mathieu Faverge's avatar
Mathieu Faverge committed
163
        chameleon_error("CHAMELEON_zhemm", "chameleon_tune() failed");
164 165 166 167
        return status;
    }

    /* Set MT & NT & KT */
168
    NB = CHAMELEON_NB;
169

Mathieu Faverge's avatar
Mathieu Faverge committed
170
    chameleon_sequence_create( chamctxt, &sequence );
171

172
    /* Submit the matrix conversion */
Mathieu Faverge's avatar
Mathieu Faverge committed
173
    chameleon_zlap2tile( chamctxt, &descAl, &descAt, ChamDescInput, uplo,
174
                     A, NB, NB, LDA, Am, Am, Am, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
175
    chameleon_zlap2tile( chamctxt, &descBl, &descBt, ChamDescInput, ChamUpperLower,
176
                     B, NB, NB, LDB, N, M,  N, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
177
    chameleon_zlap2tile( chamctxt, &descCl, &descCt, ChamDescInout, ChamUpperLower,
178
                     C, NB, NB, LDC, N, M,  N, sequence, &request );
179 180

    /* Call the tile interface */
181
    CHAMELEON_zhemm_Tile_Async(  side, uplo, alpha, &descAt, &descBt, beta, &descCt, sequence, &request );
182

Mathieu Faverge's avatar
Mathieu Faverge committed
183
    /* Submit the matrix conversion back */
Mathieu Faverge's avatar
Mathieu Faverge committed
184
    chameleon_ztile2lap( chamctxt, &descAl, &descAt,
185
                     ChamDescInput, uplo, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
186
    chameleon_ztile2lap( chamctxt, &descBl, &descBt,
187
                     ChamDescInput, ChamUpperLower, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
188
    chameleon_ztile2lap( chamctxt, &descCl, &descCt,
189
                     ChamDescInout, ChamUpperLower, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
190

Mathieu Faverge's avatar
Mathieu Faverge committed
191
    chameleon_sequence_wait( chamctxt, sequence );
Mathieu Faverge's avatar
Mathieu Faverge committed
192

Mathieu Faverge's avatar
Mathieu Faverge committed
193
    /* Cleanup the temporary data */
Mathieu Faverge's avatar
Mathieu Faverge committed
194 195 196
    chameleon_ztile2lap_cleanup( chamctxt, &descAl, &descAt );
    chameleon_ztile2lap_cleanup( chamctxt, &descBl, &descBt );
    chameleon_ztile2lap_cleanup( chamctxt, &descCl, &descCt );
197 198

    status = sequence->status;
Mathieu Faverge's avatar
Mathieu Faverge committed
199
    chameleon_sequence_destroy( chamctxt, sequence );
200 201 202
    return status;
}

203 204
/**
 ********************************************************************************
205
 *
206
 * @ingroup CHAMELEON_Complex64_t_Tile
207
 *
208 209
 *  CHAMELEON_zhemm_Tile - Performs Hermitian matrix multiplication.
 *  Tile equivalent of CHAMELEON_zhemm().
210 211 212 213 214 215 216 217 218
 *  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:
219 220
 *          = ChamLeft:      \f[ C = \alpha \times A \times B + \beta \times C \f]
 *          = ChamRight:     \f[ C = \alpha \times B \times A + \beta \times C \f]
221 222 223 224
 *
 * @param[in] uplo
 *          Specifies whether the upper or lower triangular part of
 *          the hermitian matrix A is to be referenced as follows:
225
 *          = ChamLower:     Only the lower triangular part of the
226
 *                             hermitian matrix A is to be referenced.
227
 *          = ChamUpper:     Only the upper triangular part of the
228 229 230 231 232 233
 *                             hermitian matrix A is to be referenced.
 *
 * @param[in] alpha
 *          Specifies the scalar alpha.
 *
 * @param[in] A
234
 *          A is a LDA-by-ka matrix, where ka is M when side = ChamLeft,
235 236 237 238 239 240 241 242 243 244 245 246 247 248 249
 *          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.
 *
 *******************************************************************************
 *
250
 * @retval CHAMELEON_SUCCESS successful exit
251 252 253
 *
 *******************************************************************************
 *
254 255 256 257 258
 * @sa CHAMELEON_zhemm
 * @sa CHAMELEON_zhemm_Tile_Async
 * @sa CHAMELEON_chemm_Tile
 * @sa CHAMELEON_dhemm_Tile
 * @sa CHAMELEON_shemm_Tile
259
 *
260
 */
261 262 263
int CHAMELEON_zhemm_Tile( cham_side_t side, cham_uplo_t uplo,
                      CHAMELEON_Complex64_t alpha, CHAM_desc_t *A, CHAM_desc_t *B,
                      CHAMELEON_Complex64_t beta,  CHAM_desc_t *C )
264
{
Mathieu Faverge's avatar
Mathieu Faverge committed
265
    CHAM_context_t *chamctxt;
266 267
    RUNTIME_sequence_t *sequence = NULL;
    RUNTIME_request_t request = RUNTIME_REQUEST_INITIALIZER;
268 269
    int status;

Mathieu Faverge's avatar
Mathieu Faverge committed
270 271 272
    chamctxt = chameleon_context_self();
    if (chamctxt == NULL) {
        chameleon_fatal_error("CHAMELEON_zhemm_Tile", "CHAMELEON not initialized");
273
        return CHAMELEON_ERR_NOT_INITIALIZED;
274
    }
Mathieu Faverge's avatar
Mathieu Faverge committed
275
    chameleon_sequence_create( chamctxt, &sequence );
276

277
    CHAMELEON_zhemm_Tile_Async(side, uplo, alpha, A, B, beta, C, sequence, &request );
278

279 280 281
    CHAMELEON_Desc_Flush( A, sequence );
    CHAMELEON_Desc_Flush( B, sequence );
    CHAMELEON_Desc_Flush( C, sequence );
Mathieu Faverge's avatar
Mathieu Faverge committed
282

Mathieu Faverge's avatar
Mathieu Faverge committed
283
    chameleon_sequence_wait( chamctxt, sequence );
284
    status = sequence->status;
Mathieu Faverge's avatar
Mathieu Faverge committed
285
    chameleon_sequence_destroy( chamctxt, sequence );
286 287 288
    return status;
}

289 290
/**
 ********************************************************************************
291
 *
292
 * @ingroup CHAMELEON_Complex64_t_Tile_Async
293
 *
294 295
 *  CHAMELEON_zhemm_Tile_Async - Performs Hermitian matrix multiplication.
 *  Non-blocking equivalent of CHAMELEON_zhemm_Tile().
296 297 298 299 300 301 302 303 304 305 306 307 308 309
 *  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).
 *
 *******************************************************************************
 *
310 311 312 313 314
 * @sa CHAMELEON_zhemm
 * @sa CHAMELEON_zhemm_Tile
 * @sa CHAMELEON_chemm_Tile_Async
 * @sa CHAMELEON_dhemm_Tile_Async
 * @sa CHAMELEON_shemm_Tile_Async
315
 *
316
 */
317 318 319 320
int CHAMELEON_zhemm_Tile_Async( cham_side_t side, cham_uplo_t uplo,
                            CHAMELEON_Complex64_t alpha, CHAM_desc_t *A, CHAM_desc_t *B,
                            CHAMELEON_Complex64_t beta,  CHAM_desc_t *C,
                            RUNTIME_sequence_t *sequence, RUNTIME_request_t *request )
321
{
Mathieu Faverge's avatar
Mathieu Faverge committed
322
    CHAM_context_t *chamctxt;
323

Mathieu Faverge's avatar
Mathieu Faverge committed
324 325 326
    chamctxt = chameleon_context_self();
    if (chamctxt == NULL) {
        chameleon_fatal_error("CHAMELEON_zhemm_Tile_Async", "CHAMELEON not initialized");
327
        return CHAMELEON_ERR_NOT_INITIALIZED;
328 329
    }
    if (sequence == NULL) {
Mathieu Faverge's avatar
Mathieu Faverge committed
330
        chameleon_fatal_error("CHAMELEON_zhemm_Tile_Async", "NULL sequence");
331
        return CHAMELEON_ERR_UNALLOCATED;
332 333
    }
    if (request == NULL) {
Mathieu Faverge's avatar
Mathieu Faverge committed
334
        chameleon_fatal_error("CHAMELEON_zhemm_Tile_Async", "NULL request");
335
        return CHAMELEON_ERR_UNALLOCATED;
336 337
    }
    /* Check sequence status */
338 339
    if (sequence->status == CHAMELEON_SUCCESS) {
        request->status = CHAMELEON_SUCCESS;
Mathieu Faverge's avatar
Mathieu Faverge committed
340 341
    }
    else {
Mathieu Faverge's avatar
Mathieu Faverge committed
342
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_SEQUENCE_FLUSHED);
Mathieu Faverge's avatar
Mathieu Faverge committed
343
    }
344 345

    /* Check descriptors for correctness */
Mathieu Faverge's avatar
Mathieu Faverge committed
346 347 348
    if (chameleon_desc_check(A) != CHAMELEON_SUCCESS) {
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "invalid first descriptor");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
349
    }
Mathieu Faverge's avatar
Mathieu Faverge committed
350 351 352
    if (chameleon_desc_check(B) != CHAMELEON_SUCCESS) {
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "invalid second descriptor");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
353
    }
Mathieu Faverge's avatar
Mathieu Faverge committed
354 355 356
    if (chameleon_desc_check(C) != CHAMELEON_SUCCESS) {
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "invalid third descriptor");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
357 358
    }
    /* Check input arguments */
359
    if ( (side != ChamLeft) && (side != ChamRight) ){
Mathieu Faverge's avatar
Mathieu Faverge committed
360 361
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "illegal value of side");
        return chameleon_request_fail(sequence, request, -1);
362
    }
363
    if ((uplo != ChamLower) && (uplo != ChamUpper)) {
Mathieu Faverge's avatar
Mathieu Faverge committed
364 365
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "illegal value of uplo");
        return chameleon_request_fail(sequence, request, -2);
366 367 368 369
    }

    /* Check matrices sizes */
    if ( (B->m != C->m) || (B->n != C->n) ) {
Mathieu Faverge's avatar
Mathieu Faverge committed
370 371
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "B and C must have the same size");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
372 373
    }
    if ( (A->m != A->n) ||
374 375
         ( (side == ChamLeft)  && (A->m != B->m ) ) ||
         ( (side == ChamRight) && (A->m != B->n ) ) ) {
Mathieu Faverge's avatar
Mathieu Faverge committed
376 377
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "Matrix A must be square of size M or N regarding side.");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
378 379 380 381
    }

    /* Check tiles sizes */
    if ( (B->mb != C->mb) || (B->nb != C->nb) ) {
Mathieu Faverge's avatar
Mathieu Faverge committed
382 383
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "B and C must have the same tile sizes");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
384 385
    }
    if ( (A->mb != A->nb) ||
386 387
         ( (side == ChamLeft)  && (A->mb != B->mb ) ) ||
         ( (side == ChamRight) && (A->mb != B->nb ) ) ) {
Mathieu Faverge's avatar
Mathieu Faverge committed
388 389
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "Matrix A must be square with square tiles wich fits the reagding tile size of B and C");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
390 391 392 393
    }

    /* Check submatrix starting point */
    /* if ( (B->i != C->i) || (B->j != C->j) ) { */
Mathieu Faverge's avatar
Mathieu Faverge committed
394 395
    /*     chameleon_error("CHAMELEON_zhemm_Tile_Async", "B and C submatrices doesn't match"); */
    /*     return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE); */
396 397
    /* } */
    /* if ( (A->i != A->j) ||  */
398 399
    /*          ( (side == ChamLeft)  && (A->i != B->i ) ) ||  */
    /*          ( (side == ChamRight) && (A->i != B->j ) ) ) { */
Mathieu Faverge's avatar
Mathieu Faverge committed
400 401
    /*     chameleon_error("CHAMELEON_zhemm_Tile_Async", "Submatrix A must start on diagnonal and match submatrices B and C."); */
    /*     return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE); */
402 403 404 405
    /* } */
    if( (A->i != 0) || (A->j != 0) ||
        (B->i != 0) || (B->j != 0) ||
        (C->i != 0) || (C->j != 0) ) {
Mathieu Faverge's avatar
Mathieu Faverge committed
406 407
        chameleon_error("CHAMELEON_zhemm_Tile_Async", "Submatrices are not supported for now");
        return chameleon_request_fail(sequence, request, CHAMELEON_ERR_ILLEGAL_VALUE);
408 409 410
    }

    /* Quick return */
411
    if ( (C->m == 0) || (C->n == 0) ||
412
         ( (alpha == (CHAMELEON_Complex64_t)0.0) && (beta == (CHAMELEON_Complex64_t)1.0) ) )
413
    {
414
        return CHAMELEON_SUCCESS;
415
    }
416

Mathieu Faverge's avatar
Mathieu Faverge committed
417
    chameleon_pzhemm( side, uplo, alpha, A, B, beta, C, sequence, request );
418

419
    return CHAMELEON_SUCCESS;
420
}