PLASMA  2.4.5
PLASMA - Parallel Linear Algebra for Scalable Multi-core Architectures
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups
core_sgetrf_reclap.c File Reference
#include <math.h>
#include <cblas.h>
#include <lapacke.h>
#include "common.h"
Include dependency graph for core_sgetrf_reclap.c:

Go to the source code of this file.

Macros

#define AMAX1BUF_SIZE   (48 << 1)

Functions

void CORE_sgetrf_reclap_init (void)
int CORE_sgetrf_reclap (int M, int N, float *A, int LDA, int *IPIV, int *info)
void QUARK_CORE_sgetrf_reclap (Quark *quark, Quark_Task_Flags *task_flags, int m, int n, int nb, float *A, int lda, int *IPIV, PLASMA_sequence *sequence, PLASMA_request *request, PLASMA_bool check_info, int iinfo, int nbthread)
void CORE_sgetrf_reclap_quark (Quark *quark)

Detailed Description

PLASMA core_blas kernel PLASMA is a software package provided by Univ. of Tennessee, Univ. of California Berkeley and Univ. of Colorado Denver

Version:
2.4.5
Author:
Hatem Ltaief
Mathieu Faverge
Piotr Luszczek
Date:
2009-11-15

s Tue Nov 22 14:35:21 2011

Definition in file core_sgetrf_reclap.c.


Macro Definition Documentation

#define AMAX1BUF_SIZE   (48 << 1)

Definition at line 99 of file core_sgetrf_reclap.c.


Function Documentation

int CORE_sgetrf_reclap ( int  M,
int  N,
float *  A,
int  LDA,
int *  IPIV,
int *  info 
)

Definition at line 307 of file core_sgetrf_reclap.c.

References coreblas_error, max, min, and PLASMA_SUCCESS.

{
int thidx = info[1];
int thcnt = min( info[2], M / N );
int minMN = min(M, N);
if( M < 0 ) {
coreblas_error(1, "illegal value of M");
return -1;
}
if( N < 0 ) {
coreblas_error(2, "illegal value of N");
return -2;
}
if( LDA < max(1, M) ) {
coreblas_error(5, "illegal value of LDA");
return -5;
}
/*
* Quick return
*/
if ( (M == 0) || (N == 0) || (thidx >= thcnt) ){
}
*info = 0;
CORE_sgetrf_reclap_rec( M, minMN, A, LDA, IPIV, info,
thidx, thcnt, 0 );
if ( N > minMN ) {
CORE_sgetrf_reclap_update(M, 0, minMN, N-minMN,
A, LDA, IPIV,
thidx, thcnt);
}
return info[0];
}

Here is the caller graph for this function:

void CORE_sgetrf_reclap_init ( void  )

Definition at line 105 of file core_sgetrf_reclap.c.

References AMAX1BUF_SIZE.

{
int i;
for (i = 0; i < AMAX1BUF_SIZE; ++i) CORE_samax1buf[i] = -1.0;
sfmin = LAPACKE_slamch_work('S');
}

Here is the caller graph for this function:

void CORE_sgetrf_reclap_quark ( Quark quark)

Definition at line 381 of file core_sgetrf_reclap.c.

References A, CORE_sgetrf_reclap(), IPIV, plasma_sequence_flush(), PLASMA_SUCCESS, QUARK_Get_RankInTask(), and quark_unpack_args_10.

{
int M;
int N;
float *A;
int LDA;
int *IPIV;
PLASMA_sequence *sequence;
PLASMA_request *request;
PLASMA_bool check_info;
int iinfo;
int info[3];
int maxthreads;
quark_unpack_args_10(quark, M, N, A, LDA, IPIV, sequence, request,
check_info, iinfo, maxthreads );
info[1] = QUARK_Get_RankInTask(quark);
info[2] = maxthreads;
CORE_sgetrf_reclap( M, N, A, LDA, IPIV, info );
if (info[1] == 0 && info[0] != PLASMA_SUCCESS && check_info)
plasma_sequence_flush(quark, sequence, request, iinfo + info[0] );
}

Here is the call graph for this function:

Here is the caller graph for this function:

void QUARK_CORE_sgetrf_reclap ( Quark quark,
Quark_Task_Flags task_flags,
int  m,
int  n,
int  nb,
float *  A,
int  lda,
int *  IPIV,
PLASMA_sequence sequence,
PLASMA_request request,
PLASMA_bool  check_info,
int  iinfo,
int  nbthread 
)

Definition at line 351 of file core_sgetrf_reclap.c.

References CORE_sgetrf_reclap_quark(), DAG_CORE_GETRF, INOUT, OUTPUT, QUARK_Insert_Task(), and VALUE.

{
sizeof(int), &m, VALUE,
sizeof(int), &n, VALUE,
sizeof(float)*nb*nb, A, INOUT,
sizeof(int), &lda, VALUE,
sizeof(int)*nb, IPIV, OUTPUT,
sizeof(PLASMA_sequence*), &sequence, VALUE,
sizeof(PLASMA_request*), &request, VALUE,
sizeof(PLASMA_bool), &check_info, VALUE,
sizeof(int), &iinfo, VALUE,
sizeof(int), &nbthread, VALUE,
0);
}

Here is the call graph for this function:

Here is the caller graph for this function: