Mentions légales du service

Skip to content
Snippets Groups Projects
Commit 4869799b authored by BOUCHERIE Raphael's avatar BOUCHERIE Raphael
Browse files

renamed pzgeqrfhqr in pzgeqrf_param, modified cmake for libhqr depencies,...

renamed pzgeqrfhqr in pzgeqrf_param, modified cmake for libhqr depencies, added compute file for param, need to change the commentary and add the test file
parent 4726f7e3
No related branches found
No related tags found
1 merge request!47Integration of hierarchicah householder reduction trees
...@@ -250,6 +250,7 @@ endif() ...@@ -250,6 +250,7 @@ endif()
######################### #########################
set(CHAMELEON_DEP "") set(CHAMELEON_DEP "")
set(CHAMELEON_USE_LIBHQR 1)
find_package(LIBHQR REQUIRED) find_package(LIBHQR REQUIRED)
list(INSERT CHAMELEON_DEP 0 ${LIBHQR_LIBRARIES}) list(INSERT CHAMELEON_DEP 0 ${LIBHQR_LIBRARIES})
include_directories(${LIBHQR_INCLUDE_DIRS}) include_directories(${LIBHQR_INCLUDE_DIRS})
......
...@@ -105,7 +105,7 @@ set(ZSRC ...@@ -105,7 +105,7 @@ set(ZSRC
pzgelqfrh.c pzgelqfrh.c
pzgeqrf.c pzgeqrf.c
pzgeqrfrh.c pzgeqrfrh.c
pzgeqrfhqr.c pzgeqrf_param.c
pzgetrf_incpiv.c pzgetrf_incpiv.c
pzgetrf_nopiv.c pzgetrf_nopiv.c
pzlacpy.c pzlacpy.c
...@@ -138,6 +138,7 @@ set(ZSRC ...@@ -138,6 +138,7 @@ set(ZSRC
zgelqf.c zgelqf.c
zgelqs.c zgelqs.c
zgeqrf.c zgeqrf.c
zgeqrf_param.c
zgeqrs.c zgeqrs.c
#zgesv.c #zgesv.c
zgesv_incpiv.c zgesv_incpiv.c
......
...@@ -9,7 +9,7 @@ ...@@ -9,7 +9,7 @@
/** /**
* *
* @file pzgeqrfhqr.c * @file pzgeqrf_param.c
* *
* MORSE auxiliary routines * MORSE auxiliary routines
* MORSE is a software package provided by Univ. of Tennessee, * MORSE is a software package provided by Univ. of Tennessee,
...@@ -30,26 +30,26 @@ ...@@ -30,26 +30,26 @@
#include "control/common.h" #include "control/common.h"
#include "libhqr.h" #include "libhqr.h"
#define A(m,n) A, (m), (n) #define A(m,n) A, (m), (n)
#define TS(m,n) TS, (m), (n) #define TS(m,n) TS, (m), (n)
#define TT(m,n) TT, (m), (n) #define TT(m,n) TT, (m), (n)
#if defined(CHAMELEON_COPY_DIAG) #if defined(CHAMELEON_COPY_DIAG)
#define DIAG(m,n) DIAG, (m), (n) #define D(m,n) D, (m), (n)
#else #else
#define DIAG(m,n) A, (m), (n) #define D(m,n) A, (m), (n)
#endif #endif
/***************************************************************************//** /**
* Parallel tile QR factorization (reduction Householder) - dynamic scheduling * Parallel tile QR factorization (reduction Householder) - dynamic scheduling
**/ */
void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT, void morse_pzgeqrf_param( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT,
MORSE_sequence_t *sequence, MORSE_request_t *request) MORSE_sequence_t *sequence, MORSE_request_t *request)
{ {
MORSE_context_t *morse; MORSE_context_t *morse;
MORSE_option_t options; MORSE_option_t options;
size_t ws_worker = 0; size_t ws_worker = 0;
size_t ws_host = 0; size_t ws_host = 0;
MORSE_desc_t *DIAG = NULL; MORSE_desc_t *D = NULL;
int k, m, n, i, j, p; int k, m, n, i, j, p;
int K, M, RD; int K, M, RD;
...@@ -90,23 +90,6 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -90,23 +90,6 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
tiles = (int*)malloc((qrtree->mt)*sizeof(int)); tiles = (int*)malloc((qrtree->mt)*sizeof(int));
memset( tiles, 0, (qrtree->mt)*sizeof(int) ); memset( tiles, 0, (qrtree->mt)*sizeof(int) );
#if defined(CHAMELEON_USE_MAGMA)
/* Worker space
*
* zgeqrt = max( A->nb * (ib+1), ib * (ib + A->nb) )
* ztsqrt = max( A->nb * (ib+1), ib * (ib + A->nb) )
*/
ws_worker = chameleon_max( ws_worker, ib * (ib + A->nb) );
/* Host space
*
* zgeqrt = ib * (A->nb+3*ib) + A->nb )
* ztsqrt = 2 * ib * (A->nb+ib) + A->nb
*/
ws_host = chameleon_max( ws_host, ib * (A->mb + 3 * ib) + A->mb );
ws_host = chameleon_max( ws_host, 2 * ib * (A->nb + ib) + A->nb );
#endif
ws_worker *= sizeof(MORSE_Complex64_t); ws_worker *= sizeof(MORSE_Complex64_t);
ws_host *= sizeof(MORSE_Complex64_t); ws_host *= sizeof(MORSE_Complex64_t);
...@@ -115,9 +98,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -115,9 +98,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
#if defined(CHAMELEON_COPY_DIAG) #if defined(CHAMELEON_COPY_DIAG)
{ {
/* necessary to avoid dependencies between tasks regarding the diag tile */ /* necessary to avoid dependencies between tasks regarding the diag tile */
int nblk = ( A->mt + BS -1 ) / BS; D = (MORSE_desc_t*)malloc(sizeof(MORSE_desc_t));
DIAG = (MORSE_desc_t*)malloc(sizeof(MORSE_desc_t)); morse_zdesc_alloc(*D, A->mb, A->nb, A->m, A->n, 0, 0, A->m, A->n, );
morse_zdesc_alloc_diag(*DIAG, A->mb, A->nb, nblk * A->mb, A->nb, 0, 0, nblk * A->mb, A->nb, A->p, A->q);
} }
#endif #endif
...@@ -138,7 +120,7 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -138,7 +120,7 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
MORSE_TASK_zgeqrt( MORSE_TASK_zgeqrt(
&options, &options,
tempmm, tempkn, ib, TS->nb, tempmm, tempkn, ib, TS->nb,
A(m, k), ldam, A( m, k), ldam,
TS(m, k), TS->mb); TS(m, k), TS->mb);
if ( k < (A->nt-1) ) { if ( k < (A->nt-1) ) {
#if defined(CHAMELEON_COPY_DIAG) #if defined(CHAMELEON_COPY_DIAG)
...@@ -146,13 +128,13 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -146,13 +128,13 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
&options, &options,
MorseLower, tempmm, A->nb, A->nb, MorseLower, tempmm, A->nb, A->nb,
A(m, k), ldam, A(m, k), ldam,
DIAG(m, k), ldam ); D(m, k), ldam );
#if defined(CHAMELEON_USE_CUDA) #if defined(CHAMELEON_USE_CUDA)
MORSE_TASK_zlaset( MORSE_TASK_zlaset(
&options, &options,
MorseUpper, tempmm, A->nb, MorseUpper, tempmm, A->nb,
0., 1., 0., 1.,
DIAG(m, k), ldam ); D(m, k), ldam );
#endif #endif
#endif #endif
} }
...@@ -162,9 +144,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -162,9 +144,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
&options, &options,
MorseLeft, MorseConjTrans, MorseLeft, MorseConjTrans,
tempmm, tempnn, tempkmin, ib, TS->nb, tempmm, tempnn, tempkmin, ib, TS->nb,
DIAG(m, k), ldam, D( m, k), ldam,
TS(m, k), TS->mb, TS(m, k), TS->mb,
A(m, n), ldam); A( m, n), ldam);
} }
} }
...@@ -184,8 +166,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -184,8 +166,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
MORSE_TASK_ztsqrt( MORSE_TASK_ztsqrt(
&options, &options,
tempmm, tempkn, ib, TS->nb, tempmm, tempkn, ib, TS->nb,
A(p, k), ldap, A( p, k), ldap,
A(m, k), ldam, A( m, k), ldam,
TS(m, k), TS->mb); TS(m, k), TS->mb);
for (n = k+1; n < A->nt; n++) { for (n = k+1; n < A->nt; n++) {
...@@ -194,9 +176,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -194,9 +176,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
&options, &options,
MorseLeft, MorseConjTrans, MorseLeft, MorseConjTrans,
A->nb, tempnn, tempmm, tempnn, A->nb, ib, TS->nb, A->nb, tempnn, tempmm, tempnn, A->nb, ib, TS->nb,
A(p, n), ldap, A( p, n), ldap,
A(m, n), ldam, A( m, n), ldam,
A(m, k), ldam, A( m, k), ldam,
TS(m, k), TS->mb); TS(m, k), TS->mb);
} }
} }
...@@ -206,8 +188,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -206,8 +188,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
MORSE_TASK_zttqrt( MORSE_TASK_zttqrt(
&options, &options,
tempmm, tempkn, ib, TT->nb, tempmm, tempkn, ib, TT->nb,
A(p, k), ldap, A( p, k), ldap,
A(m, k), ldam, A( m, k), ldam,
TT(m, k), TT->mb); TT(m, k), TT->mb);
for (n = k+1; n < A->nt; n++) { for (n = k+1; n < A->nt; n++) {
...@@ -216,9 +198,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -216,9 +198,9 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
&options, &options,
MorseLeft, MorseConjTrans, MorseLeft, MorseConjTrans,
A->mb, tempnn, tempmm, tempnn, A->nb, ib, TT->nb, A->mb, tempnn, tempmm, tempnn, A->nb, ib, TT->nb,
A(p, n), ldap, A( p, n), ldap,
A(m, n), ldam, A( m, n), ldam,
A(m, k), ldam, A( m, k), ldam,
TT(m, k), TT->mb); TT(m, k), TT->mb);
} }
} }
...@@ -233,8 +215,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_ ...@@ -233,8 +215,8 @@ void morse_pzgeqrfhqr( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_
#if defined(CHAMELEON_COPY_DIAG) #if defined(CHAMELEON_COPY_DIAG)
MORSE_Sequence_Wait(sequence); MORSE_Sequence_Wait(sequence);
morse_desc_mat_free(DIAG); morse_desc_mat_free(D);
free(DIAG); free(D);
#endif #endif
(void)DIAG; (void)D;
} }
...@@ -205,7 +205,8 @@ int MORSE_zgeqrf_Tile(MORSE_desc_t *A, MORSE_desc_t *T) ...@@ -205,7 +205,8 @@ int MORSE_zgeqrf_Tile(MORSE_desc_t *A, MORSE_desc_t *T)
return status; return status;
} }
/***************************************************************************//** /**
*****************************************************************************
* *
* @ingroup MORSE_Complex64_t_Tile_Async * @ingroup MORSE_Complex64_t_Tile_Async
* *
......
/**
*
* @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 zgeqrf_param.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 Jakub Kurzak
* @author Mathieu Faverge
* @author Emmanuel Agullo
* @author Cedric Castagnede
* @date 2010-11-15
* @precisions normal z -> s d c
*
**/
#include "control/common.h"
/**
*******************************************************************************
*
* @ingroup MORSE_Complex64_t
*
* MORSE_zgeqrf_param - Computes the tile QR factorization of a complex M-by-N
* matrix A: A = Q * R.
*
*******************************************************************************
*
* @param[in] M
* The number of rows of the matrix A. M >= 0.
*
* @param[in] N
* The number of columns of the matrix A. N >= 0.
*
* @param[in] A
* On entry, the M-by-N matrix A.
* On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N
* upper trapezoidal matrix R (R is upper triangular if M >= N); the elements below the
* diagonal represent the unitary matrix Q as a product of elementary reflectors stored
* by tiles.
*
* @param[in] LDA
* The leading dimension of the array A. LDA >= max(1,M).
*
* @param[in] descTS
* On exit, auxiliary factorization data, required by MORSE_zgeqrs to solve the system
* of equations.
*
*******************************************************************************
*
* @return
* \retval MORSE_SUCCESS successful exit
* \retval <0 if -i, the i-th argument had an illegal value
*
*******************************************************************************
*
* @sa MORSE_zgeqrf_param_Tile
* @sa MORSE_zgeqrf_param_Tile_Async
* @sa MORSE_cgeqrf
* @sa MORSE_dgeqrf
* @sa MORSE_sgeqrf
* @sa MORSE_zgeqrs
*
******************************************************************************/
int MORSE_zgeqrf_param(const libhqr_tree_t *qrtree, int M, int N,
MORSE_Complex64_t *A, int LDA,
MORSE_desc_t *descTS, MORSE_desc_t *descTT)
{
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_zgeqrf_param", "MORSE not initialized");
return MORSE_ERR_NOT_INITIALIZED;
}
/* Check input arguments */
if (M < 0) {
morse_error("MORSE_zgeqrf_param", "illegal value of M");
return -1;
}
if (N < 0) {
morse_error("MORSE_zgeqrf_param", "illegal value of N");
return -2;
}
if (LDA < chameleon_max(1, M)) {
morse_error("MORSE_zgeqrf_param", "illegal value of LDA");
return -4;
}
/* Quick return */
if (chameleon_min(M, N) == 0)
return MORSE_SUCCESS;
/* Tune NB & IB depending on M, N & NRHS; Set NBNBSIZE */
status = morse_tune(MORSE_FUNC_ZGELS, M, N, 0);
if (status != MORSE_SUCCESS) {
morse_error("MORSE_zgeqrf_param", "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, M, N, sequence, &request,
morse_desc_mat_free(&(descA)) );
/* } else {*/
/* morse_ziplap2tile( descA, A, NB, NB, LDA, N, 0, 0, M, N,*/
/* sequence, &request);*/
/* }*/
/* Call the tile interface */
MORSE_zgeqrf_param_Tile_Async(qrtree, &descA, descTS, descTT, 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_zgeqrf_param_Tile - Computes the tile QR factorization of a matrix.
* Tile equivalent of MORSE_zgeqrf_param().
* Operates on matrices stored by tiles.
* All matrices are passed through descriptors.
* All dimensions are taken from the descriptors.
*
*******************************************************************************
*
* @param[in,out] A
* On entry, the M-by-N matrix A.
* On exit, the elements on and above the diagonal of the array contain the min(M,N)-by-N
* upper trapezoidal matrix R (R is upper triangular if M >= N); the elements below the
* diagonal represent the unitary matrix Q as a product of elementary reflectors stored
* by tiles.
*
* @param[out] T
* On exit, auxiliary factorization data, required by MORSE_zgeqrs to solve the system
* of equations.
*
*******************************************************************************
*
* @return
* \retval MORSE_SUCCESS successful exit
*
*******************************************************************************
*
* @sa MORSE_zgeqrf_param
* @sa MORSE_zgeqrf_param_Tile_Async
* @sa MORSE_cgeqrf_param_Tile
* @sa MORSE_dgeqrf_param_Tile
* @sa MORSE_sgeqrf_param_Tile
* @sa MORSE_zgeqrs_param_Tile
*
******************************************************************************/
int MORSE_zgeqrf_param_Tile(const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT)
{
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_zgeqrf_param_Tile", "MORSE not initialized");
return MORSE_ERR_NOT_INITIALIZED;
}
morse_sequence_create(morse, &sequence);
MORSE_zgeqrf_param_Tile_Async(qrtree, A, TS, TT, sequence, &request);
morse_sequence_wait(morse, sequence);
RUNTIME_desc_getoncpu(A);
status = sequence->status;
morse_sequence_destroy(morse, sequence);
return status;
}
/**
*****************************************************************************
*
* @ingroup MORSE_Complex64_t_Tile_Async
*
* MORSE_zgeqrf_param_Tile_Async - Computes the tile QR factorization of a matrix.
* Non-blocking equivalent of MORSE_zgeqrf_param_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_zgeqrf_param
* @sa MORSE_zgeqrf_param_Tile
* @sa MORSE_cgeqrf_param_Tile_Async
* @sa MORSE_dgeqrf_param_Tile_Async
* @sa MORSE_sgeqrf_param_Tile_Async
* @sa MORSE_zgeqrs_param_Tile_Async
*
******************************************************************************/
int MORSE_zgeqrf_param_Tile_Async(const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT,
MORSE_sequence_t *sequence, MORSE_request_t *request)
{
MORSE_context_t *morse;
morse = morse_context_self();
if (morse == NULL) {
morse_error("MORSE_zgeqrf_param_Tile", "MORSE not initialized");
return MORSE_ERR_NOT_INITIALIZED;
}
if (sequence == NULL) {
morse_fatal_error("MORSE_zgeqrf_param_Tile", "NULL sequence");
return MORSE_ERR_UNALLOCATED;
}
if (request == NULL) {
morse_fatal_error("MORSE_zgeqrf_param_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);
/* Check descriptors for correctness */
if (morse_desc_check(A) != MORSE_SUCCESS) {
morse_error("MORSE_zgeqrf_param_Tile", "invalid first descriptor");
return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
}
if (morse_desc_check(TS) != MORSE_SUCCESS) {
morse_error("MORSE_zgeqrf_param_Tile", "invalid second descriptor");
return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
}
if (morse_desc_check(TT) != MORSE_SUCCESS) {
morse_error("MORSE_zgeqrf_param_Tile", "invalid second descriptor");
return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
}
/* Check input arguments */
if (A->nb != A->mb) {
morse_error("MORSE_zgeqrf_param_Tile", "only square tiles supported");
return morse_request_fail(sequence, request, MORSE_ERR_ILLEGAL_VALUE);
}
/* Quick return */
/*
if (chameleon_min(M, N) == 0)
return MORSE_SUCCESS;
*/
morse_pzgeqrf_param(qrtree, A, TS, TT, sequence, request);
return MORSE_SUCCESS;
}
...@@ -154,3 +154,8 @@ void morse_pzunmqrrh(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_d ...@@ -154,3 +154,8 @@ void morse_pzunmqrrh(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_d
void morse_pzunmlq(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_desc_t *B, MORSE_desc_t *T, MORSE_sequence_t *sequence, MORSE_request_t *request); void morse_pzunmlq(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_desc_t *B, MORSE_desc_t *T, MORSE_sequence_t *sequence, MORSE_request_t *request);
void morse_pzunmlqrh(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_desc_t *B, MORSE_desc_t *T, int BS, MORSE_sequence_t *sequence, MORSE_request_t *request); void morse_pzunmlqrh(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, MORSE_desc_t *B, MORSE_desc_t *T, int BS, MORSE_sequence_t *sequence, MORSE_request_t *request);
void morse_pzbuild( MORSE_enum uplo, MORSE_desc_t *A, void *user_data, void* user_build_callback, MORSE_sequence_t *sequence, MORSE_request_t *request ); void morse_pzbuild( MORSE_enum uplo, MORSE_desc_t *A, void *user_data, void* user_build_callback, MORSE_sequence_t *sequence, MORSE_request_t *request );
#if defined(CHAMELEON_USE_LIBHQR)
void morse_pzgeqrf_param( const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT,
MORSE_sequence_t *sequence, MORSE_request_t *request);
#endif /* defined(CHAMELEON_USE_LIBHQR) */
...@@ -37,6 +37,9 @@ ...@@ -37,6 +37,9 @@
/* Communication engine */ /* Communication engine */
#cmakedefine CHAMELEON_USE_MPI #cmakedefine CHAMELEON_USE_MPI
/* Householder reduction trees for QR like operations */
#cmakedefine CHAMELEON_USE_LIBHQR
/* GPU Support */ /* GPU Support */
#cmakedefine CHAMELEON_USE_CUDA #cmakedefine CHAMELEON_USE_CUDA
#cmakedefine CHAMELEON_USE_CUBLAS #cmakedefine CHAMELEON_USE_CUBLAS
...@@ -45,5 +48,4 @@ ...@@ -45,5 +48,4 @@
/* Simulating */ /* Simulating */
#cmakedefine CHAMELEON_SIMULATION #cmakedefine CHAMELEON_SIMULATION
#endif /* CHAMELEON_CONFIG_H_HAS_BEEN_INCLUDED */ #endif /* CHAMELEON_CONFIG_H_HAS_BEEN_INCLUDED */
...@@ -32,7 +32,6 @@ ...@@ -32,7 +32,6 @@
#define MORSE_VERSION_MINOR @MORSE_VERSION_MINOR@ #define MORSE_VERSION_MINOR @MORSE_VERSION_MINOR@
#define MORSE_VERSION_MICRO @MORSE_VERSION_MICRO@ #define MORSE_VERSION_MICRO @MORSE_VERSION_MICRO@
/* **************************************************************************** /* ****************************************************************************
* MORSE types and constants * MORSE types and constants
*/ */
...@@ -122,6 +121,10 @@ int MORSE_Sequence_Wait (MORSE_sequence_t *sequence); ...@@ -122,6 +121,10 @@ int MORSE_Sequence_Wait (MORSE_sequence_t *sequence);
} }
#endif #endif
#if defined(CHAMELEON_USE_LIBHQR)
#include "libhqr.h"
#endif /* defined(CHAMELEON_USE_LIBHQR) */
#include "morse_z.h" #include "morse_z.h"
#include "morse_c.h" #include "morse_c.h"
#include "morse_d.h" #include "morse_d.h"
......
...@@ -273,9 +273,17 @@ int MORSE_zunmqr_Tile_Async(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A, ...@@ -273,9 +273,17 @@ int MORSE_zunmqr_Tile_Async(MORSE_enum side, MORSE_enum trans, MORSE_desc_t *A,
//int MORSE_zgecfi_Async(int m, int n, MORSE_Complex64_t *A, MORSE_enum f_in, int imb, int inb, MORSE_enum f_out, int omb, int onb, MORSE_sequence_t *sequence, MORSE_request_t *request); //int MORSE_zgecfi_Async(int m, int n, MORSE_Complex64_t *A, MORSE_enum f_in, int imb, int inb, MORSE_enum f_out, int omb, int onb, MORSE_sequence_t *sequence, MORSE_request_t *request);
//int MORSE_zgetmi_Async(int m, int n, MORSE_Complex64_t *A, MORSE_enum f_in, int mb, int inb, MORSE_sequence_t *sequence, MORSE_request_t *request); //int MORSE_zgetmi_Async(int m, int n, MORSE_Complex64_t *A, MORSE_enum f_in, int mb, int inb, MORSE_sequence_t *sequence, MORSE_request_t *request);
/** **************************************************************************** /**
* Declarations of libhqr dependent functions.
*/
#if defined(CHAMELEON_USE_LIBHQR)
int MORSE_zgeqrf_param(const libhqr_tree_t *qrtree, int M, int N, MORSE_Complex64_t *A, int LDA, MORSE_desc_t *descTS, MORSE_desc_t *descTT);
int MORSE_zgeqrf_param_Tile(const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT);
int MORSE_zgeqrf_param_Tile_Async(const libhqr_tree_t *qrtree, MORSE_desc_t *A, MORSE_desc_t *TS, MORSE_desc_t *TT, MORSE_sequence_t *sequence, MORSE_request_t *request);
#endif /* defined(CHAMELEON_USE_LIBHQR) */
/**
* Declarations of workspace allocation functions (tile layout) - alphabetical order * Declarations of workspace allocation functions (tile layout) - alphabetical order
**/ */
int MORSE_Alloc_Workspace_zgesv_incpiv( int N, MORSE_desc_t **descL, int **IPIV, int p, int q); int MORSE_Alloc_Workspace_zgesv_incpiv( int N, MORSE_desc_t **descL, int **IPIV, int p, int q);
int MORSE_Alloc_Workspace_zgetrf_incpiv(int M, int N, MORSE_desc_t **descL, int **IPIV, int p, int q); int MORSE_Alloc_Workspace_zgetrf_incpiv(int M, int N, MORSE_desc_t **descL, int **IPIV, int p, int q);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment