18 #include "common/common.h"
27 #if defined( PASTIX_WITH_MPI )
52 pastix_complex64_t *b;
58 assert( colnbr <= solvmtx->colmax );
59 assert( cblk->
cblktype & CBLK_RECV );
61 #if defined (PASTIX_DEBUG_MPI)
62 fprintf( stderr,
"[%2d] RHS Bwd: Post Isend cblk %ld to %2d at index %ld of size %ld\n",
64 (
long)cblk->
lcolidx, (
long)colnbr );
72 assert( rhsb->
cblkb[ idx ] == NULL );
74 rhsb->
cblkb[ idx ] = malloc( size *
sizeof(pastix_complex64_t) );
76 rc = LAPACKE_zlacpy_work( LAPACK_COL_MAJOR,
'A', colnbr, rhsb->
n,
77 b, rhsb->
ld, rhsb->
cblkb[ idx ], colnbr );
80 b = rhsb->
cblkb[ idx ];
83 rc = MPI_Isend( b, size, PASTIX_MPI_COMPLEX64,
85 assert( rc == MPI_SUCCESS );
87 solverCommMatrixAdd( solvmtx, cblk->
ownerid, size *
sizeof(pastix_complex64_t) );
90 pastix_atomic_lock( &(solvmtx->
reqlock) );
93 assert( solvmtx->
reqnum >= 0 );
100 pastix_atomic_unlock( &(solvmtx->
reqlock) );
127 cpucblk_zrequest_rhs_bwd_handle_send(
const args_solve_t *enums,
134 assert( cblk->
cblktype & CBLK_RECV );
135 assert( enums->solve_step == PastixSolveBackward );
137 #if defined(PASTIX_DEBUG_MPI)
139 size_t cblksize =
cblk_colnbr( cblk ) * rhsb->
n *
sizeof(pastix_complex64_t);
141 fprintf( stderr,
"[%2d] RHS Bwd: Isend for cblk %ld toward %2d ( %ld Bytes ) (DONE)\n",
142 solvmtx->clustnum, (
long)cblk->
gcblknum, cblk->
ownerid, (
long)cblksize );
146 if ( rhsb->
cblkb[ idx ] ) {
147 memFree_null( rhsb->
cblkb[ idx ] );
181 cpucblk_zrequest_rhs_bwd_handle_recv(
const args_solve_t *enums,
185 const MPI_Status *status,
186 pastix_complex64_t *recvbuf )
189 int src = status->MPI_SOURCE;
190 int tag = status->MPI_TAG;
193 assert( ( 0 <= src ) && ( src < solvmtx->clustnbr ) );
194 assert( ( 0 <= tag ) && ( tag < solvmtx->gcblknbr ) );
200 assert( cblk->
cblktype & CBLK_FANIN );
202 #if defined(PASTIX_DEBUG_MPI)
210 size = colnbr * rhsb->
n *
sizeof(pastix_complex64_t);
212 rc = MPI_Get_count( status, MPI_CHAR, &count );
213 assert( rc == MPI_SUCCESS );
214 assert( count == size );
217 fprintf( stderr,
"[%2d] RHS Bwd : recv of size %d/%ld for cblk %ld (DONE)\n",
218 solvmtx->clustnum, count, (
long)size, (
long)cblk->
gcblknum );
223 rhsb->
cblkb[ idx ] = recvbuf;
226 if ( solvmtx->computeQueue ) {
236 assert( rhsb->
cblkb[ idx ] == NULL );
278 cpucblk_zrequest_rhs_bwd_handle(
const args_solve_t *enums,
284 const MPI_Status *statuses )
287 int nbrequest = outcount;
289 for ( i = 0; i < outcount; i++ ) {
295 if ( solvmtx->
reqidx[reqid] == -1 ) {
297 pastix_complex64_t *recvbuf;
301 memcpy( &status, statuses + i,
sizeof(MPI_Status) );
302 MPI_Get_count( &status, PASTIX_MPI_COMPLEX64, &size );
304 MALLOC_INTERN( recvbuf, size, pastix_complex64_t );
305 memcpy( recvbuf, solvmtx->
rcoeftab, size *
sizeof(pastix_complex64_t) );
312 MPI_Start( solvmtx->
reqtab + reqid );
316 MPI_Request_free( solvmtx->
reqtab + reqid );
317 solvmtx->
reqtab[reqid] = MPI_REQUEST_NULL;
320 cpucblk_zrequest_rhs_bwd_handle_recv( enums, solvmtx, rhsb,
321 threadid, &status, recvbuf );
328 assert( cblk->
cblktype & CBLK_RECV );
330 cpucblk_zrequest_rhs_bwd_handle_send( enums, solvmtx, rhsb, cblk );
333 solvmtx->
reqidx[ reqid ] = -1;
369 cpucblk_zmpi_rhs_bwd_progress(
const args_solve_t *enums,
374 pthread_t tid = pthread_self();
376 int nbrequest, nbfree;
377 int indexes[ solvmtx->
reqnbr ];
378 MPI_Status statuses[ solvmtx->
reqnbr ];
381 pthread_mutex_lock( &pastix_comm_lock );
382 if ( pastix_comm_tid == (pthread_t)-1 ) {
383 pastix_comm_tid = tid;
385 pthread_mutex_unlock( &pastix_comm_lock );
387 if ( tid != pastix_comm_tid ) {
401 pastix_atomic_lock( &(solvmtx->
reqlock) );
402 nbrequest = solvmtx->
reqnum;
403 pastix_atomic_unlock( &(solvmtx->
reqlock) );
405 while( (outcount > 0) && (nbrequest > 0) )
407 MPI_Testsome( nbrequest, solvmtx->
reqtab, &outcount, indexes, statuses );
411 if ( outcount > 0 ) {
412 nbfree = cpucblk_zrequest_rhs_bwd_handle( enums, solvmtx, rhsb, threadid,
413 outcount, indexes, statuses );
420 pastix_atomic_lock( &(solvmtx->
reqlock) );
422 cpucblk_zupdate_reqtab( solvmtx );
424 nbrequest = solvmtx->
reqnum;
425 pastix_atomic_unlock( &(solvmtx->
reqlock) );
428 pastix_comm_tid = (pthread_t)-1;
469 #if defined(PASTIX_WITH_MPI)
470 if ( cblk->
cblktype & CBLK_FANIN ) {
485 cpucblk_zmpi_rhs_bwd_progress( enums, solvmtx, rhsb, rank );
488 assert( !(cblk->
cblktype & (CBLK_FANIN | CBLK_RECV)) );
489 do { pastix_yield(); }
while( cblk->
ctrbcnt > 0 );
532 ctrbcnt = pastix_atomic_dec_32b( &(fcbk->
ctrbcnt) );
534 #if defined(PASTIX_WITH_MPI)
536 cpucblk_zisend_rhs_bwd( solvmtx, rhsb, fcbk );
542 if ( solvmtx->computeQueue ) {
591 #if defined(PASTIX_WITH_MPI)
595 int reqnbr = solvmtx->
reqnum;
598 #if defined(PASTIX_DEBUG_MPI)
599 fprintf( stderr,
"[%2d] Wait for all pending communications\n",
603 for ( i=0; i<reqnbr; i++ ) {
604 if ( solvmtx->
reqtab[i] == MPI_REQUEST_NULL ) {
611 assert( solvmtx->
reqidx[i] != -1 );
613 rc = MPI_Wait( solvmtx->
reqtab + i, &status );
614 assert( rc == MPI_SUCCESS );
619 assert( cblk->
cblktype & CBLK_RECV );
621 cpucblk_zrequest_rhs_bwd_handle_send( enums, solvmtx, rhsb, cblk );
625 assert( solvmtx->
reqnum == 0 );
BEGIN_C_DECLS typedef int pastix_int_t
static void pqueuePush1(pastix_queue_t *q, pastix_int_t elt, double key1)
Push an element with a single key.
void cpucblk_zrelease_rhs_bwd_deps(const args_solve_t *enums, SolverMatrix *solvmtx, pastix_rhs_t rhsb, const SolverCblk *cblk, SolverCblk *fcbk)
Release the dependencies of the given cblk after an update.
void cpucblk_zrequest_rhs_bwd_cleanup(const args_solve_t *enums, pastix_int_t sched, SolverMatrix *solvmtx, pastix_rhs_t rhsb)
Waitall routine for current cblk request.
int cpucblk_zincoming_rhs_bwd_deps(int rank, const args_solve_t *enums, SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Wait for incoming dependencies, and return when cblk->ctrbcnt has reached 0.
void solve_cblk_ztrsmsp_backward(const args_solve_t *enums, SolverMatrix *datacode, SolverCblk *cblk, pastix_rhs_t rhsb)
Apply a backward solve related to one cblk to all the right hand side.
Main PaStiX RHS structure.
pastix_atomic_lock_t reqlock
static pastix_int_t cblk_colnbr(const SolverCblk *cblk)
Compute the number of columns in a column block.
SolverBlok *restrict bloktab
pastix_int_t *restrict browtab
SolverCblk *restrict cblktab
Solver column block structure.
Solver column block structure.