Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
solverstack
PaStiX
Commits
f2d6647c
Commit
f2d6647c
authored
Jul 07, 2015
by
Mathieu Faverge
Browse files
Add infrastructure for muti-threaded code
parent
f726dd4e
Changes
10
Hide whitespace changes
Inline
Side-by-side
pastix/cmake_modules/morse/precision_generator/subs.py
View file @
f2d6647c
...
...
@@ -145,6 +145,8 @@ subs = {
(
''
,
'csc_s'
,
'csc_d'
,
'csc_c'
,
'csc_z'
),
(
''
,
'sequential_s'
,
'sequential_d'
,
'sequential_c'
,
'sequential_z'
),
(
''
,
'coeftab_s'
,
'coeftab_d'
,
'coeftab_c'
,
'coeftab_z'
),
(
''
,
'thread_s'
,
'thread_d'
,
'thread_c'
,
'thread_z'
),
(
''
,
'thread_ps'
,
'thread_pd'
,
'thread_pc'
,
'thread_pz'
),
# ----- Complex numbers
# \b regexp here avoids conjugate -> conjfugate,
...
...
pastix/sopalin/pastix_task_solve.c
View file @
f2d6647c
...
...
@@ -26,27 +26,27 @@
#include "s_bcsc.h"
void
pastix_
static_
trsm
(
int
flttype
,
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
void
*
b
,
int
ldb
)
pastix_trsm
(
int
flttype
,
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
void
*
b
,
int
ldb
)
{
switch
(
flttype
)
{
case
PastixComplex64
:
pastix_static
_ztrsm
(
side
,
uplo
,
trans
,
diag
,
sequential
_ztrsm
(
side
,
uplo
,
trans
,
diag
,
sopalin_data
,
nrhs
,
(
pastix_complex64_t
*
)
b
,
ldb
);
break
;
case
PastixComplex32
:
pastix_static
_ctrsm
(
side
,
uplo
,
trans
,
diag
,
sequential
_ctrsm
(
side
,
uplo
,
trans
,
diag
,
sopalin_data
,
nrhs
,
(
pastix_complex32_t
*
)
b
,
ldb
);
break
;
case
PastixDouble
:
trans
=
(
trans
==
PastixConjTrans
)
?
PastixTrans
:
trans
;
pastix_static
_dtrsm
(
side
,
uplo
,
trans
,
diag
,
sequential
_dtrsm
(
side
,
uplo
,
trans
,
diag
,
sopalin_data
,
nrhs
,
(
double
*
)
b
,
ldb
);
break
;
case
PastixFloat
:
trans
=
(
trans
==
PastixConjTrans
)
?
PastixTrans
:
trans
;
pastix_static
_strsm
(
side
,
uplo
,
trans
,
diag
,
sequential
_strsm
(
side
,
uplo
,
trans
,
diag
,
sopalin_data
,
nrhs
,
(
float
*
)
b
,
ldb
);
break
;
default:
...
...
@@ -55,9 +55,9 @@ pastix_static_trsm( int flttype, int side, int uplo, int trans, int diag,
}
void
pastix_
static_
diag
(
int
flttype
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
void
*
b
,
int
ldb
)
pastix_diag
(
int
flttype
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
void
*
b
,
int
ldb
)
{
switch
(
flttype
)
{
case
PastixComplex64
:
...
...
@@ -176,26 +176,26 @@ pastix_task_solve( pastix_data_t *pastix_data,
clockStart
(
timer
);
switch
(
pastix_data
->
iparm
[
IPARM_FACTORIZATION
]
){
case
PastixFactLLT
:
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixConjTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixConjTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
break
;
case
PastixFactLDLT
:
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
diag
(
pastix_data
->
bcsc
->
flttype
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_diag
(
pastix_data
->
bcsc
->
flttype
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
break
;
case
PastixFactLDLH
:
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
diag
(
pastix_data
->
bcsc
->
flttype
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixConjTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_diag
(
pastix_data
->
bcsc
->
flttype
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixConjTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
break
;
case
PastixFactLU
:
default:
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_
static_
trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixUpper
,
PastixNoTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixLower
,
PastixNoTrans
,
PastixUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
pastix_trsm
(
pastix_data
->
bcsc
->
flttype
,
PastixLeft
,
PastixUpper
,
PastixNoTrans
,
PastixNonUnit
,
&
sopalin_data
,
nrhs
,
b
,
ldb
);
break
;
}
clockStop
(
timer
);
...
...
pastix/sopalin/pastix_task_sopalin.c
View file @
f2d6647c
...
...
@@ -15,36 +15,17 @@
*
**/
#include "common.h"
#include "isched.h"
#include "csc.h"
#include "bcsc.h"
#include "sopalin_data.h"
static
void
(
*
sopalinFacto
[
4
][
4
])(
sopalin_data_t
*
)
=
{
{
pastix_static_spotrf
,
pastix_static_dpotrf
,
pastix_static_cpotrf
,
pastix_static_zpotrf
},
{
pastix_static_ssytrf
,
pastix_static_dsytrf
,
pastix_static_csytrf
,
pastix_static_zsytrf
},
{
pastix_static_sgetrf
,
pastix_static_dgetrf
,
pastix_static_cgetrf
,
pastix_static_zgetrf
},
{
pastix_static_ssytrf
,
pastix_static_dsytrf
,
pastix_static_chetrf
,
pastix_static_zhetrf
}
{
sopalin_spotrf
,
sopalin_dpotrf
,
sopalin_cpotrf
,
sopalin_zpotrf
},
{
sopalin_ssytrf
,
sopalin_dsytrf
,
sopalin_csytrf
,
sopalin_zsytrf
},
{
sopalin_sgetrf
,
sopalin_dgetrf
,
sopalin_cgetrf
,
sopalin_zgetrf
},
{
sopalin_ssytrf
,
sopalin_dsytrf
,
sopalin_chetrf
,
sopalin_zhetrf
}
};
void
...
...
@@ -172,6 +153,7 @@ int
pastix_task_sopalin
(
pastix_data_t
*
pastix_data
,
pastix_csc_t
*
csc
)
{
extern
isched_t
*
scheduler
;
sopalin_data_t
sopalin_data
;
SolverBackup_t
*
sbackup
;
/* #ifdef PASTIX_WITH_MPI */
...
...
@@ -255,6 +237,7 @@ pastix_task_sopalin( pastix_data_t *pastix_data,
sopalin_data
.
diagthreshold
=
pastix_data
->
dparm
[
DPARM_EPSILON_MAGN_CTRL
]
*
pastix_data
->
dparm
[
DPARM_A_NORM
];
}
}
sopalin_data
.
sched
=
scheduler
;
sbackup
=
solverBackupInit
(
pastix_data
->
solvmatr
);
pastix_data
->
solvmatr
->
restore
=
2
;
...
...
pastix/sopalin/sequential_zdiag.c
View file @
f2d6647c
...
...
@@ -16,6 +16,7 @@
* @precisions normal z -> s d c
*
**/
#include <cblas.h>
#include "common.h"
#include "sopalin_data.h"
#include "pastix_zcores.h"
...
...
pastix/sopalin/sequential_zgetrf.c
View file @
f2d6647c
/**
*
* @file s
equential
_zgetrf.c
* @file s
opalin
_zgetrf.c
*
* PaStiX factorization routines
* PaStiX is a software package provided by Inria Bordeaux - Sud-Ouest,
...
...
@@ -17,11 +17,12 @@
*
**/
#include "common.h"
#include "isched.h"
#include "sopalin_data.h"
#include "pastix_zcores.h"
void
pastix_static
_zgetrf
(
sopalin_data_t
*
sopalin_data
)
sequential
_zgetrf
(
sopalin_data_t
*
sopalin_data
)
{
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
...
...
@@ -37,3 +38,60 @@ pastix_static_zgetrf( sopalin_data_t *sopalin_data )
coeftab_zdump
(
datacode
,
"getrf_L.txt"
);
#endif
}
void
thread_pzgetrf
(
int
rank
,
void
*
args
)
{
sopalin_data_t
*
sopalin_data
=
(
sopalin_data_t
*
)
args
;
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
Task
*
t
;
pastix_int_t
i
,
ii
;
pastix_int_t
tasknbr
,
*
tasktab
;
tasknbr
=
datacode
->
ttsknbr
[
rank
];
tasktab
=
datacode
->
ttsktab
[
rank
];
for
(
ii
=
0
;
ii
<
tasknbr
;
ii
++
)
{
i
=
tasktab
[
ii
];
t
=
datacode
->
tasktab
+
i
;
cblk
=
datacode
->
cblktab
+
t
->
cblknum
;
/* Compute */
core_zgetrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
}
#if defined(PASTIX_DEBUG_FACTO)
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
if
(
rank
==
0
)
{
coeftab_zdump
(
datacode
,
"getrf_L.txt"
);
}
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
#endif
}
void
thread_zgetrf
(
sopalin_data_t
*
sopalin_data
)
{
isched_parallel_call
(
sopalin_data
->
sched
,
thread_pzgetrf
,
sopalin_data
);
}
static
void
(
*
zgetrf_table
[
4
])(
sopalin_data_t
*
)
=
{
sequential_zgetrf
,
thread_zgetrf
,
NULL
,
NULL
};
void
sopalin_zgetrf
(
sopalin_data_t
*
sopalin_data
)
{
int
sched
=
0
;
void
(
*
zgetrf
)(
sopalin_data_t
*
)
=
zgetrf_table
[
sched
];
if
(
zgetrf
==
NULL
)
{
zgetrf
=
thread_zgetrf
;
}
zgetrf
(
sopalin_data
);
}
pastix/sopalin/sequential_zhetrf.c
View file @
f2d6647c
/**
*
* @file s
equential
_zhetrf.c
* @file s
opalin
_zhetrf.c
*
* PaStiX factorization routines
* PaStiX is a software package provided by Inria Bordeaux - Sud-Ouest,
...
...
@@ -17,11 +17,12 @@
*
**/
#include "common.h"
#include "isched.h"
#include "sopalin_data.h"
#include "pastix_zcores.h"
void
pastix_static
_zhetrf
(
sopalin_data_t
*
sopalin_data
)
sequential
_zhetrf
(
sopalin_data_t
*
sopalin_data
)
{
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
...
...
@@ -37,3 +38,61 @@ pastix_static_zhetrf( sopalin_data_t *sopalin_data )
coeftab_zdump
(
datacode
,
"hetrf_L.txt"
);
#endif
}
void
thread_pzhetrf
(
int
rank
,
void
*
args
)
{
sopalin_data_t
*
sopalin_data
=
(
sopalin_data_t
*
)
args
;
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
Task
*
t
;
pastix_int_t
i
,
ii
;
pastix_int_t
tasknbr
,
*
tasktab
;
tasknbr
=
datacode
->
ttsknbr
[
rank
];
tasktab
=
datacode
->
ttsktab
[
rank
];
for
(
ii
=
0
;
ii
<
tasknbr
;
ii
++
)
{
i
=
tasktab
[
ii
];
t
=
datacode
->
tasktab
+
i
;
cblk
=
datacode
->
cblktab
+
t
->
cblknum
;
/* Compute */
core_zhetrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
}
#if defined(PASTIX_DEBUG_FACTO)
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
if
(
rank
==
0
)
{
coeftab_zdump
(
datacode
,
"hetrf_L.txt"
);
}
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
#endif
}
void
thread_zhetrf
(
sopalin_data_t
*
sopalin_data
)
{
isched_parallel_call
(
sopalin_data
->
sched
,
thread_pzhetrf
,
sopalin_data
);
}
static
void
(
*
zhetrf_table
[
4
])(
sopalin_data_t
*
)
=
{
sequential_zhetrf
,
thread_zhetrf
,
NULL
,
NULL
};
void
sopalin_zhetrf
(
sopalin_data_t
*
sopalin_data
)
{
int
sched
=
0
;
void
(
*
zhetrf
)(
sopalin_data_t
*
)
=
zhetrf_table
[
sched
];
if
(
zhetrf
==
NULL
)
{
zhetrf
=
thread_zhetrf
;
}
zhetrf
(
sopalin_data
);
}
pastix/sopalin/sequential_zpotrf.c
View file @
f2d6647c
/**
*
* @file s
equential
_zpotrf.c
* @file s
opalin
_zpotrf.c
*
* PaStiX factorization routines
* PaStiX is a software package provided by Inria Bordeaux - Sud-Ouest,
...
...
@@ -17,11 +17,12 @@
*
**/
#include "common.h"
#include "isched.h"
#include "sopalin_data.h"
#include "pastix_zcores.h"
void
pastix_static
_zpotrf
(
sopalin_data_t
*
sopalin_data
)
sequential
_zpotrf
(
sopalin_data_t
*
sopalin_data
)
{
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
...
...
@@ -30,10 +31,67 @@ pastix_static_zpotrf( sopalin_data_t *sopalin_data )
cblk
=
datacode
->
cblktab
;
for
(
i
=
0
;
i
<
datacode
->
cblknbr
;
i
++
,
cblk
++
){
/* Compute */
core_zpotrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
core_zpotrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
}
#if defined(PASTIX_DEBUG_FACTO)
coeftab_zdump
(
datacode
,
"potrf_L.txt"
);
#endif
}
void
thread_pzpotrf
(
int
rank
,
void
*
args
)
{
sopalin_data_t
*
sopalin_data
=
(
sopalin_data_t
*
)
args
;
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
Task
*
t
;
pastix_int_t
i
,
ii
;
pastix_int_t
tasknbr
,
*
tasktab
;
tasknbr
=
datacode
->
ttsknbr
[
rank
];
tasktab
=
datacode
->
ttsktab
[
rank
];
for
(
ii
=
0
;
ii
<
tasknbr
;
ii
++
)
{
i
=
tasktab
[
ii
];
t
=
datacode
->
tasktab
+
i
;
cblk
=
datacode
->
cblktab
+
t
->
cblknum
;
/* Compute */
core_zpotrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
}
#if defined(PASTIX_DEBUG_FACTO)
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
if
(
rank
==
0
)
{
coeftab_zdump
(
datacode
,
"potrf_L.txt"
);
}
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
#endif
}
void
thread_zpotrf
(
sopalin_data_t
*
sopalin_data
)
{
isched_parallel_call
(
sopalin_data
->
sched
,
thread_pzpotrf
,
sopalin_data
);
}
static
void
(
*
zpotrf_table
[
4
])(
sopalin_data_t
*
)
=
{
sequential_zpotrf
,
thread_zpotrf
,
NULL
,
NULL
};
void
sopalin_zpotrf
(
sopalin_data_t
*
sopalin_data
)
{
int
sched
=
0
;
void
(
*
zpotrf
)(
sopalin_data_t
*
)
=
zpotrf_table
[
sched
];
if
(
zpotrf
==
NULL
)
{
zpotrf
=
thread_zpotrf
;
}
zpotrf
(
sopalin_data
);
}
pastix/sopalin/sequential_zsytrf.c
View file @
f2d6647c
/**
*
* @file s
equential
_zsytrf.c
* @file s
opalin
_zsytrf.c
*
* PaStiX factorization routines
* PaStiX is a software package provided by Inria Bordeaux - Sud-Ouest,
...
...
@@ -17,11 +17,12 @@
*
**/
#include "common.h"
#include "isched.h"
#include "sopalin_data.h"
#include "pastix_zcores.h"
void
pastix_static
_zsytrf
(
sopalin_data_t
*
sopalin_data
)
sequential
_zsytrf
(
sopalin_data_t
*
sopalin_data
)
{
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
...
...
@@ -37,3 +38,60 @@ pastix_static_zsytrf( sopalin_data_t *sopalin_data )
coeftab_zdump
(
datacode
,
"sytrf_L.txt"
);
#endif
}
void
thread_pzsytrf
(
int
rank
,
void
*
args
)
{
sopalin_data_t
*
sopalin_data
=
(
sopalin_data_t
*
)
args
;
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
;
Task
*
t
;
pastix_int_t
i
,
ii
;
pastix_int_t
tasknbr
,
*
tasktab
;
tasknbr
=
datacode
->
ttsknbr
[
rank
];
tasktab
=
datacode
->
ttsktab
[
rank
];
for
(
ii
=
0
;
ii
<
tasknbr
;
ii
++
)
{
i
=
tasktab
[
ii
];
t
=
datacode
->
tasktab
+
i
;
cblk
=
datacode
->
cblktab
+
t
->
cblknum
;
/* Compute */
core_zsytrfsp1d
(
datacode
,
cblk
,
sopalin_data
->
diagthreshold
);
}
#if defined(PASTIX_DEBUG_FACTO)
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
if
(
rank
==
0
)
{
coeftab_zdump
(
datacode
,
"sytrf_L.txt"
);
}
isched_barrier_wait
(
&
(((
isched_t
*
)(
sopalin_data
->
sched
))
->
barrier
)
);
#endif
}
void
thread_zsytrf
(
sopalin_data_t
*
sopalin_data
)
{
isched_parallel_call
(
sopalin_data
->
sched
,
thread_pzsytrf
,
sopalin_data
);
}
static
void
(
*
zsytrf_table
[
4
])(
sopalin_data_t
*
)
=
{
sequential_zsytrf
,
thread_zsytrf
,
NULL
,
NULL
};
void
sopalin_zsytrf
(
sopalin_data_t
*
sopalin_data
)
{
int
sched
=
0
;
void
(
*
zsytrf
)(
sopalin_data_t
*
)
=
zsytrf_table
[
sched
];
if
(
zsytrf
==
NULL
)
{
zsytrf
=
thread_zsytrf
;
}
zsytrf
(
sopalin_data
);
}
pastix/sopalin/sequential_ztrsm.c
View file @
f2d6647c
...
...
@@ -26,9 +26,9 @@ static pastix_complex64_t zone = 1.0;
static
pastix_complex64_t
mzone
=
-
1
.
0
;
void
pastix_static
_ztrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex64_t
*
b
,
int
ldb
)
sequential
_ztrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex64_t
*
b
,
int
ldb
)
{
SolverMatrix
*
datacode
=
sopalin_data
->
solvmtx
;
SolverCblk
*
cblk
,
*
fcbk
;
...
...
pastix/sopalin/sopalin_data.h
View file @
f2d6647c
...
...
@@ -18,38 +18,44 @@
#define _SOPALIN_DATA_H_
struct
sopalin_data_s
{
void
*
sched
;
SolverMatrix
*
solvmtx
;
double
diagthreshold
;
/* Threshold for static pivoting on diagonal value */
};
typedef
struct
sopalin_data_s
sopalin_data_t
;
void
pastix_static_ztrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex64_t
*
b
,
int
ldb
);
void
pastix_static_ctrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex32_t
*
b
,
int
ldb
);
void
pastix_static_dtrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
double
*
b
,
int
ldb
);
void
pastix_static_strsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
float
*
b
,
int
ldb
);
void
coeftab_zdump
(
const
SolverMatrix
*
solvmtx
,
const
char
*
filename
);
void
coeftab_cdump
(
const
SolverMatrix
*
solvmtx
,
const
char
*
filename
);
void
coeftab_ddump
(
const
SolverMatrix
*
solvmtx
,
const
char
*
filename
);
void
coeftab_sdump
(
const
SolverMatrix
*
solvmtx
,
const
char
*
filename
);
void
sequential_ztrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex64_t
*
b
,
int
ldb
);
void
sequential_ctrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex32_t
*
b
,
int
ldb
);
void
sequential_dtrsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
double
*
b
,
int
ldb
);
void
sequential_strsm
(
int
side
,
int
uplo
,
int
trans
,
int
diag
,
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
float
*
b
,
int
ldb
);
void
sequential_zdiag
(
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex64_t
*
b
,
int
ldb
);
void
sequential_cdiag
(
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
pastix_complex32_t
*
b
,
int
ldb
);
void
sequential_ddiag
(
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
double
*
b
,
int
ldb
);
void
sequential_sdiag
(
sopalin_data_t
*
sopalin_data
,
int
nrhs
,
float
*
b
,
int
ldb
);
void
pastix_static
_zgetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_cgetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_dgetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_sgetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_zgetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_cgetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_dgetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_sgetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_zhetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_chetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_zhetrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_chetrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_zpotrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_cpotrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_dpotrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_spotrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_zpotrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_cpotrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_dpotrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_spotrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_zsytrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_csytrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_dsytrf
(
sopalin_data_t
*
sopalin_data
);
void
pastix_static
_ssytrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_zsytrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_csytrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_dsytrf
(
sopalin_data_t
*
sopalin_data
);
void
sopalin
_ssytrf
(
sopalin_data_t
*
sopalin_data
);
#endif
/* _SOPALIN_DATA_H_ */
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment