zbuild.c 9.71 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
/**
 *
 * @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 zbuild.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
 * @author Guillaume Sylvand
 * @date 2016-09-05
 * @precisions normal z -> s d c
 *
 **/
#include "control/common.h"

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
 *
 * @ingroup MORSE_Complex64_t
 *
 *  MORSE_zbuild - Generate a matrix by calling user provided function.
 *
 *******************************************************************************
 *
 * @param[in] uplo
 *          Specifies the part of the matrix A to be copied to B.
 *            = MorseUpperLower: All the matrix A
 *            = MorseUpper: Upper triangular part
 *            = MorseLower: Lower triangular part
 *
 * @param[in] M
 *          The number of rows of A.
 *
 * @param[in] N
 *          The order of the matrix A. N >= 0.
 *
 * @param[out] A
 *          On exit, The matrix A generated.
 *
 * @param[in] LDA
 *          The leading dimension of the array A. LDA >= max(1,M).
 *
 * @param[in] user_data
 *          The user data used in the matrix generation, it will be passed by chameleon
 *          to the build callback function (see below).
 *
 * @param[in] user_build_callback
 *          The user function to call to generate tiles.
 *          The prototype of the callback is :
 *          void myFcn(int row_min, int row_max, int col_min, int col_max, void *buffer, int ld, void *user_data)
 *          It is expected to build the block of matrix [row_min, row_max] x [col_min, col_max]
 *          (with both min and max values included in the intervals,
 *          index start at 0 like in C, NOT 1 like in Fortran)
 *          and store it at the adresse 'buffer' with leading dimension 'ld'
 *          The argument 'user_data' is an opaque pointer on any user data, it is passed by
 *          the user to Morse_zbuild (see above) and transmitted by chameleon to the callback.
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *          \retval <0 if -i, the i-th argument had an illegal value
 *
 *******************************************************************************
 *
 * @sa MORSE_zbuild_Tile
 * @sa MORSE_zbuild_Tile_Async
 * @sa MORSE_cbuild
 * @sa MORSE_dbuild
 * @sa MORSE_sbuild
 *
 ******************************************************************************/
int MORSE_zbuild( MORSE_enum uplo, int M, int N,
                  MORSE_Complex64_t *A, int LDA,
                  void *user_data, void* user_build_callback)
{
Mathieu Faverge's avatar
Mathieu Faverge committed
89 90 91 92 93 94
    int NB;
    int status;
    MORSE_context_t *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t request = MORSE_REQUEST_INITIALIZER;
    MORSE_desc_t descA;
95

Mathieu Faverge's avatar
Mathieu Faverge committed
96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116
    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zbuild", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    /* Check input arguments */
    if (M < 0) {
        morse_error("MORSE_zbuild", "illegal value of M");
        return -1;
    }
    if (N < 0) {
        morse_error("MORSE_zbuild", "illegal value of N");
        return -2;
    }
    if (LDA < chameleon_max(1, M)) {
        morse_error("MORSE_zbuild", "illegal value of LDA");
        return -4;
    }
    /* Quick return */
    if (chameleon_min(M, N) == 0)
        return MORSE_SUCCESS;
117

Mathieu Faverge's avatar
Mathieu Faverge committed
118 119 120 121 122 123
    /* Tune NB depending on M, N & NRHS; Set NBNB */
    status = morse_tune(MORSE_FUNC_ZGEMM, M, N, 0);
    if (status != MORSE_SUCCESS) {
        morse_error("MORSE_zbuild", "morse_tune() failed");
        return status;
    }
124

Mathieu Faverge's avatar
Mathieu Faverge committed
125 126 127 128
    /* Set NT */
    NB = MORSE_NB;
    morse_sequence_create(morse, &sequence);
    morse_zdesc_alloc(descA, NB, NB, LDA, N, 0, 0, M, N, morse_desc_mat_free(&descA));
129

Mathieu Faverge's avatar
Mathieu Faverge committed
130 131
    /* Call the tile interface */
    MORSE_zbuild_Tile_Async(uplo, &descA, user_data, user_build_callback, sequence, &request );
132

Mathieu Faverge's avatar
Mathieu Faverge committed
133 134 135
    morse_ztile2lap( morse, &descAl, &descAt,
                     MorseUpperLower, sequence, &request );

Mathieu Faverge's avatar
Mathieu Faverge committed
136
    morse_sequence_wait(morse, sequence);
Mathieu Faverge's avatar
Mathieu Faverge committed
137

Mathieu Faverge's avatar
Mathieu Faverge committed
138
    morse_desc_mat_free(&descA);
139

Mathieu Faverge's avatar
Mathieu Faverge committed
140 141
    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
142

Mathieu Faverge's avatar
Mathieu Faverge committed
143
    return status;
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
 *
 * @ingroup MORSE_Complex64_t_Tile
 *
 *  MORSE_zbuild_Tile - Generate a matrix by tiles by calling user provided function.
 *  Tile equivalent of MORSE_zbuild().
 *  Operates on matrices stored by tiles.
 *  All matrices are passed through descriptors.
 *  All dimensions are taken from the descriptors.
 *
 *******************************************************************************
 *
 * @param[in] uplo
 *          Specifies the part of the matrix A to be copied to B.
 *            = MorseUpperLower: All the matrix A
 *            = MorseUpper: Upper triangular part
 *            = MorseLower: Lower triangular part
 *
 * @param[in] A
 *          On exit, The matrix A generated.
 *
 * @param[in] user_data
 *          The data used in the matrix generation.
 *
 * @param[in] user_build_callback
 *          The function called by the codelet to fill the tiles (see MORSE_zbuild)
 *
 *******************************************************************************
 *
 * @return
 *          \retval MORSE_SUCCESS successful exit
 *
 *******************************************************************************
 *
 * @sa MORSE_zbuild
 * @sa MORSE_zbuild_Tile_Async
 * @sa MORSE_cbuild_Tile
 * @sa MORSE_dbuild_Tile
 * @sa MORSE_sbuild_Tile
 *
 ******************************************************************************/
int MORSE_zbuild_Tile( MORSE_enum uplo, MORSE_desc_t *A,
                       void *user_data, void* user_build_callback )
{
Mathieu Faverge's avatar
Mathieu Faverge committed
191 192 193 194
    MORSE_context_t *morse;
    MORSE_sequence_t *sequence = NULL;
    MORSE_request_t request = MORSE_REQUEST_INITIALIZER;
    int status;
195

Mathieu Faverge's avatar
Mathieu Faverge committed
196 197 198 199 200 201 202
    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zbuild_Tile", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    morse_sequence_create(morse, &sequence);
    MORSE_zbuild_Tile_Async( uplo, A, user_data, user_build_callback, sequence, &request );
Mathieu Faverge's avatar
Mathieu Faverge committed
203

Mathieu Faverge's avatar
Mathieu Faverge committed
204
    morse_sequence_wait(morse, sequence);
Mathieu Faverge's avatar
Mathieu Faverge committed
205

Mathieu Faverge's avatar
Mathieu Faverge committed
206 207 208
    status = sequence->status;
    morse_sequence_destroy(morse, sequence);
    return status;
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
 *
 * @ingroup MORSE_Complex64_t_Tile_Async
 *
 *  MORSE_zbuild_Tile_Async - Generate a matrix by tiles by calling user provided function.
 *  Non-blocking equivalent of MORSE_zbuild_Tile().
 *  May return before the computation is finished.
 *  Allows for pipelining of operations at runtime.
 *
 *******************************************************************************
 *
 * @param[in] uplo
 *          Specifies the part of the matrix A to be copied to B.
 *            = MorseUpperLower: All the matrix A
 *            = MorseUpper: Upper triangular part
 *            = MorseLower: Lower triangular part
 *
 * @param[in] A
 *          On exit, The matrix A generated.
 *
 * @param[in] user_data
 *          The data used in the matrix generation.
 *
 * @param[in] user_build_callback
 *          The function called by the codelet to fill the tiles (see MORSE_zbuild)
 *
 * @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_zbuild
 * @sa MORSE_zbuild_Tile
 * @sa MORSE_cbuild_Tile_Async
 * @sa MORSE_dbuild_Tile_Async
 * @sa MORSE_sbuild_Tile_Async
 *
 ******************************************************************************/
int MORSE_zbuild_Tile_Async( MORSE_enum uplo, MORSE_desc_t     *A,
                             void *user_data, void* user_build_callback,
                             MORSE_sequence_t *sequence,
                             MORSE_request_t  *request)
{
Mathieu Faverge's avatar
Mathieu Faverge committed
259
    MORSE_context_t *morse;
260

Mathieu Faverge's avatar
Mathieu Faverge committed
261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278
    morse = morse_context_self();
    if (morse == NULL) {
        morse_fatal_error("MORSE_zbuild_Tile", "MORSE not initialized");
        return MORSE_ERR_NOT_INITIALIZED;
    }
    if (sequence == NULL) {
        morse_fatal_error("MORSE_zbuild_Tile", "NULL sequence");
        return MORSE_ERR_UNALLOCATED;
    }
    if (request == NULL) {
        morse_fatal_error("MORSE_zbuild_Tile", "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);
279

Mathieu Faverge's avatar
Mathieu Faverge committed
280 281 282 283 284
    /* Check descriptors for correctness */
    if (morse_desc_check(A) != MORSE_SUCCESS) {
        morse_error("MORSE_zbuild_Tile", "invalid descriptor");
        return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
    }
285

Mathieu Faverge's avatar
Mathieu Faverge committed
286 287 288
    /* Quick return */
    if (chameleon_min( A->m, A->n ) == 0)
        return MORSE_SUCCESS;
289

Mathieu Faverge's avatar
Mathieu Faverge committed
290
    morse_pzbuild(uplo, A, user_data, user_build_callback, sequence,  request);
291

Mathieu Faverge's avatar
Mathieu Faverge committed
292
    return MORSE_SUCCESS;
293
}