PaStiX Handbook  6.3.1
cpucblk_dmpi_rhs_fwd.c
Go to the documentation of this file.
1 /**
2  *
3  * @file cpucblk_dmpi_rhs_fwd.c
4  *
5  * Precision dependent routines to manag communications for the solve part.
6  *
7  * @copyright 2015-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8  * Univ. Bordeaux. All rights reserved.
9  *
10  * @version 6.3.1
11  * @author Pierre Ramet
12  * @author Mathieu Faverge
13  * @author Tony Delarue
14  * @author Alycia Lisito
15  * @date 2023-09-20
16  *
17  * @generated from /builds/solverstack/pastix/kernels/cpucblk_zmpi_rhs_fwd.c, normal z -> d, Thu Nov 23 09:55:48 2023
18  *
19  **/
20 #include "common/common.h"
21 #include "common/pastixdata.h"
22 #include <lapacke.h>
23 #include "blend/solver.h"
25 #include "kernels.h"
26 #include "pastix_dcores.h"
27 #include "pastix_dlrcores.h"
28 
29 /**
30  *******************************************************************************
31  *
32  * @brief Send the rhs associated to a cblk->lcolidx to the remote node.
33  *
34  *******************************************************************************
35  *
36  * @param[in] solvmtx
37  * The solver matrix holding the communicator.
38  *
39  * @param[in] cblk
40  * The cblk which defines the part to sent.
41  *
42  * @param[in] rhsb
43  * The pointer to the rhs data structure that holds the vectors of the
44  * right hand side.
45  *
46  *******************************************************************************/
47 void
49  SolverCblk *cblk,
50  pastix_rhs_t rhsb )
51 {
52 #if defined(PASTIX_WITH_MPI)
53  double *b;
54  pastix_int_t colnbr = cblk_colnbr(cblk);
55  pastix_int_t size = colnbr * rhsb->n;
56  pastix_int_t idx = - cblk->bcscnum - 1;
57  int rc;
58 
59  assert( colnbr <= solvmtx->colmax );
60  assert( cblk->cblktype & CBLK_FANIN );
61  assert( rhsb->cblkb[ idx ] != NULL );
62 
63 #if defined (PASTIX_DEBUG_MPI)
64  fprintf( stderr, "[%2d] RHS Fwd: Send cblk %ld to %2d at index %ld of size %ld\n",
65  solvmtx->clustnum, (long)cblk->gcblknum, cblk->ownerid,
66  (long)cblk->lcolidx, (long)colnbr );
67 #endif
68 
69  b = (double*)(rhsb->cblkb[ idx ]);
70  assert( b != NULL );
71 
72  rc = MPI_Send( b, size, PASTIX_MPI_DOUBLE,
73  cblk->ownerid, cblk->gcblknum, solvmtx->solv_comm );
74  assert( rc == MPI_SUCCESS );
75 
76  memFree_null( rhsb->cblkb[ idx ] );
77 
78  (void)rc;
79 #else
80  (void)solvmtx;
81  (void)cblk;
82  (void)rhsb;
83 #endif
84 }
85 
86 /**
87  *******************************************************************************
88  *
89  * @brief Send the rhs associated to a cblk->lcolidx to the remote node.
90  *
91  *******************************************************************************
92  *
93  * @param[in] solvmtx
94  * The solver matrix holding the communicator.
95  *
96  * @param[in] cblk
97  * The cblk which defines the part to sent.
98  *
99  * @param[in] rhsb
100  * The pointer to the rhs data structure that holds the vectors of the
101  * right hand side.
102  *
103  *******************************************************************************/
104 void
106  SolverCblk *cblk,
107  pastix_rhs_t rhsb )
108 {
109 #if defined(PASTIX_WITH_MPI)
110  double *b = rhsb->b;
111  pastix_int_t colnbr = cblk_colnbr(cblk);
112  pastix_int_t idx = - cblk->bcscnum - 1;
113  int rc;
114 
115  assert( colnbr <= solvmtx->colmax );
116  assert( cblk->cblktype & CBLK_RECV );
117  assert( rhsb->cblkb[ idx ] == NULL );
118 
119 #if defined (PASTIX_DEBUG_MPI)
120  fprintf( stderr, "[%2d] RHS Bwd: Send cblk %ld to %2d at index %ld of size %ld\n",
121  solvmtx->clustnum, (long)cblk->gcblknum, cblk->ownerid,
122  (long)cblk->lcolidx, (long)colnbr );
123 #endif
124 
125  b += cblk->lcolidx;
126  if ( rhsb->n > 1 ) {
127  rhsb->cblkb[ idx ] = malloc( colnbr * rhsb->n * sizeof(double) );
128  rc = LAPACKE_dlacpy_work( LAPACK_COL_MAJOR, 'A', colnbr, rhsb->n, b, rhsb->ld, rhsb->cblkb[ idx ], colnbr );
129  assert( rc == 0 );
130  b = rhsb->cblkb[ idx ];
131  }
132 
133  rc = MPI_Send( b, colnbr * rhsb->n, PASTIX_MPI_DOUBLE,
134  cblk->ownerid, cblk->gcblknum, solvmtx->solv_comm );
135  assert( rc == MPI_SUCCESS );
136 
137  if ( rhsb->n > 1 ) {
138  memFree_null( rhsb->cblkb[ idx ] );
139  }
140 
141  (void)rc;
142 #else
143  (void)solvmtx;
144  (void)cblk;
145  (void)rhsb;
146 #endif
147 }
148 
149 /**
150  *******************************************************************************
151  *
152  * @brief Receive the rhs associated to a cblk->lcolidx to the remote node.
153  *
154  *******************************************************************************
155  *
156  * @param[in] solvmtx
157  * The solver matrix holding the communicator.
158  *
159  * @param[in] cblk
160  * The cblk which may define the part to sent.
161  *
162  * @param[inout] rhsb
163  * The pointer to the rhs data structure that holds the vectors of the
164  * right hand side.
165  *
166  *******************************************************************************/
167 void
169  SolverCblk *cblk,
170  pastix_rhs_t rhsb )
171 {
172 #if defined(PASTIX_WITH_MPI)
173  MPI_Status status;
174  pastix_int_t colnbr = cblk_colnbr(cblk);
175  pastix_int_t idx = - cblk->bcscnum - 1;
176  int rc;
177 
178  assert( colnbr <= solvmtx->colmax );
179  assert( cblk->cblktype & CBLK_FANIN );
180 
181 #if defined (PASTIX_DEBUG_MPI)
182  fprintf( stderr, "[%2d] RHS Bwd: Recv cblk %ld from %ld at index %ld of size %ld\n",
183  solvmtx->clustnum, (long)cblk->gcblknum, (long)cblk->ownerid,
184  (long)cblk->lcolidx, (long)colnbr );
185 #endif
186 
187  assert( rhsb->cblkb[ idx ] == NULL );
188  rhsb->cblkb[ idx ] = malloc( colnbr * rhsb->n * sizeof(double) );
189 
190  rc = MPI_Recv( rhsb->cblkb[ idx ], colnbr * rhsb->n, PASTIX_MPI_DOUBLE,
191  cblk->ownerid, cblk->gcblknum, solvmtx->solv_comm, &status );
192  assert( rc == MPI_SUCCESS );
193 
194 #if defined (PASTIX_DEBUG_MPI)
195  fprintf( stderr, "[%2d] RHS Bwd: Received cblk %ld from %2d\n",
196  solvmtx->clustnum, (long)cblk->gcblknum, status.MPI_SOURCE );
197 #endif
198 
199  (void)rc;
200 #else
201  (void)solvmtx;
202  (void)cblk;
203  (void)rhsb;
204 #endif
205 }
206 
207 /**
208  *******************************************************************************
209  *
210  * @brief Receive the rhs associated to a cblk->lcolidx to the remote node.
211  *
212  *******************************************************************************
213  *
214  * @param[in] solvmtx
215  * The solver matrix holding the communicator.
216  *
217  * @param[in] cblk
218  * The cblk which may define the part to sent.
219  *
220  * @param[inout] work
221  * The temporary buffer to receive the remote data
222  *
223  * @param[inout] rhsb
224  * The pointer to the rhs data structure that holds the vectors of the
225  * right hand side.
226  *
227  *******************************************************************************/
228 void
230  SolverCblk *cblk,
231  double *work,
232  pastix_rhs_t rhsb )
233 {
234 #if defined(PASTIX_WITH_MPI)
235  double *b = rhsb->b;
236  pastix_int_t colnbr = cblk_colnbr(cblk);
237  MPI_Status status;
238  int rc;
239 
240  assert( colnbr <= solvmtx->colmax );
241  assert( cblk->cblktype & CBLK_RECV );
242 
243 #if defined (PASTIX_DEBUG_MPI)
244  fprintf( stderr, "[%2d] RHS Fwd: Recv cblk %ld from %ld at index %ld of size %ld\n",
245  solvmtx->clustnum, (long)cblk->gcblknum, (long)cblk->ownerid,
246  (long)cblk->lcolidx, (long)colnbr );
247 #endif
248 
249  rc = MPI_Recv( work, colnbr * rhsb->n, PASTIX_MPI_DOUBLE,
250  cblk->ownerid, cblk->gcblknum, solvmtx->solv_comm, &status );
251  assert( rc == MPI_SUCCESS );
252 
253 #if defined (PASTIX_DEBUG_MPI)
254  fprintf( stderr, "[%2d] RHS Fwd: Received cblk %ld from %2d\n",
255  solvmtx->clustnum, (long)cblk->gcblknum, status.MPI_SOURCE );
256 #endif
257 
258  b += cblk->lcolidx;
259  core_dgeadd( PastixNoTrans, colnbr, rhsb->n,
260  1., work, colnbr,
261  1., b, rhsb->ld );
262 
263  (void)rc;
264 #else
265  (void)solvmtx;
266  (void)cblk;
267  (void)work;
268  (void)rhsb;
269 #endif
270 }
271 
272 #if defined( PASTIX_WITH_MPI )
273 /**
274  *******************************************************************************
275  *
276  * @brief Asynchronously sends the rhs associated to cblk->ownerid
277  *
278  *******************************************************************************
279  *
280  * @param[inout] solvmtx
281  * The solver matrix structure.
282  *
283  * @param[in] rhsb
284  * The pointer to the rhs data structure that holds the vectors of the
285  * right hand side.
286  *
287  * @param[in] cblk
288  * The column block that will be sent.
289  *
290  *******************************************************************************/
291 static inline void
292 cpucblk_disend_rhs_fwd( SolverMatrix *solvmtx,
293  pastix_rhs_t rhsb,
294  SolverCblk *cblk )
295 {
296  MPI_Request request;
297  double *b;
298  pastix_int_t colnbr = cblk_colnbr(cblk);
299  pastix_int_t size = colnbr * rhsb->n;
300  pastix_int_t idx = - cblk->bcscnum - 1;
301  int rc;
302 
303  assert( colnbr <= solvmtx->colmax );
304  assert( cblk->cblktype & CBLK_FANIN );
305  assert( rhsb->cblkb[ idx ] != NULL );
306 
307 #if defined(PASTIX_DEBUG_MPI)
308  fprintf( stderr, "[%2d] RHS Fwd: Post Isend cblk %ld to %2d at index %ld of size %ld\n",
309  solvmtx->clustnum, (long)cblk->gcblknum, cblk->ownerid,
310  (long)cblk->lcolidx, (long)(size * sizeof(double) ) );
311 #endif
312 
313  b = (double*)(rhsb->cblkb[ idx ]);
314  assert( b != NULL );
315 
316  rc = MPI_Isend( b, size, PASTIX_MPI_DOUBLE, cblk->ownerid, cblk->gcblknum,
317  solvmtx->solv_comm, &request );
318  assert( rc == MPI_SUCCESS );
319 
320  solverCommMatrixAdd( solvmtx, cblk->ownerid, size * sizeof(double) );
321 
322  /* Register the request to make it progress */
323  pastix_atomic_lock( &(solvmtx->reqlock) );
324 
325  assert( solvmtx->reqidx[ solvmtx->reqnum ] == -1 );
326  assert( solvmtx->reqnum >= 0 );
327  assert( solvmtx->reqnum < solvmtx->reqnbr );
328 
329  solvmtx->reqtab[ solvmtx->reqnum ] = request;
330  solvmtx->reqidx[ solvmtx->reqnum ] = cblk - solvmtx->cblktab;
331  solvmtx->reqnum++;
332 
333  pastix_atomic_unlock( &(solvmtx->reqlock) );
334 
335  (void)rc;
336 }
337 
338 /**
339  *******************************************************************************
340  *
341  * @brief Handle a finished request on a fanin
342  *
343  *******************************************************************************
344  *
345  * @param[in] enums
346  * Enums needed for the solve.
347  *
348  * @param[in] solvmtx
349  * The solver matrix structure.
350  *
351  * @param[in] rhsb
352  * The pointer to the rhs data structure that holds the vectors of the
353  * right hand side.
354  *
355  * @param[inout] cblk
356  * The cblk concerned by the computation.
357  *
358  *******************************************************************************/
359 void
360 cpucblk_drequest_rhs_fwd_handle_send( const args_solve_t *enums,
361  SolverMatrix *solvmtx,
362  pastix_rhs_t rhsb,
363  const SolverCblk *cblk )
364 {
365  pastix_int_t idx = - cblk->bcscnum - 1;
366 
367  assert( cblk->cblktype & CBLK_FANIN );
368  assert( enums->solve_step == PastixSolveForward );
369 
370 #if defined(PASTIX_DEBUG_MPI)
371  {
372  size_t cblksize = cblk_colnbr( cblk ) * rhsb->n * sizeof(double);
373 
374  fprintf( stderr, "[%2d] RHS Fwd: Isend for cblk %ld toward %2d ( %ld Bytes ) (DONE)\n",
375  solvmtx->clustnum, (long)cblk->gcblknum, cblk->ownerid, (long)cblksize );
376  }
377 #endif
378 
379  memFree_null( rhsb->cblkb[ idx ] );
380  (void)solvmtx;
381 }
382 
383 /**
384  *******************************************************************************
385  *
386  * @brief Handle a finished request on a recv cblk.
387  *
388  *******************************************************************************
389  *
390  * @param[in] enums
391  * Enums needed for the solve.
392  *
393  * @param[inout] solvmtx
394  * The solver matrix structure.
395  *
396  * @param[in] rhsb
397  * The pointer to the rhs data structure that holds the vectors of the
398  * right hand side.
399  *
400  * @param[inout] threadid
401  * Id of the thread calling this method.
402  *
403  * @param[inout] status
404  * Statuses of the completed request.
405  *
406  * @param[inout] recvbuf
407  * Reception buffer of the completed request.
408  *
409  *******************************************************************************/
410 static inline void
411 cpucblk_drequest_rhs_fwd_handle_recv( const args_solve_t *enums,
412  SolverMatrix *solvmtx,
413  pastix_rhs_t rhsb,
414  int threadid,
415  const MPI_Status *status,
416  double *recvbuf )
417 {
418  SolverCblk *cblk, *fcbk;
419  pastix_int_t colnbr;
420  int src = status->MPI_SOURCE;
421  int tag = status->MPI_TAG;
422  int nrhs = rhsb->n;
423 
424  assert( ( 0 <= src ) && ( src < solvmtx->clustnbr ) );
425  assert( ( 0 <= tag ) && ( tag < solvmtx->gcblknbr ) );
426 
427  /*
428  * Let's look for the local cblk
429  */
430  fcbk = solvmtx->cblktab + solvmtx->gcbl2loc[ tag ];
431  cblk = fcbk-1;
432 
433  /* Get through source */
434  while( cblk->ownerid != src ) {
435  cblk--;
436  assert( cblk >= solvmtx->cblktab );
437  assert( cblk->gcblknum == tag );
438  assert( cblk->cblktype & CBLK_RECV );
439  }
440  assert( fcbk == (solvmtx->cblktab + cblk->fblokptr->fcblknm) );
441 
442  colnbr = cblk_colnbr( cblk );
443 #if defined(PASTIX_DEBUG_MPI)
444  {
445  int rc;
446  int count = 0;
447  pastix_int_t size = colnbr * nrhs * sizeof(double);
448 
449  rc = MPI_Get_count( status, MPI_CHAR, &count );
450  assert( rc == MPI_SUCCESS );
451  assert( count == size );
452 
453  /* We can't know the sender easily, so we don't print it */
454  fprintf( stderr, "[%2d] RHS Fwd : recv of size %d/%ld for cblk %ld (DONE)\n",
455  solvmtx->clustnum, count, (long)size, (long)cblk->gcblknum );
456  }
457 #endif
458 
459  /* Initialize the cblk with the reception buffer */
460  cblk->threadid = (fcbk->threadid == -1) ? threadid : fcbk->threadid;
461 
462  {
463  double *b = rhsb->b;
464  b += cblk->lcolidx;
465  pastix_cblk_lock( fcbk );
466  core_dgeadd( PastixNoTrans, colnbr, nrhs,
467  1., recvbuf, colnbr,
468  1., b, rhsb->ld );
469  pastix_cblk_unlock( fcbk );
470  }
471 
472  /* Receptions cblks contribute to themselves */
473  cpucblk_drelease_rhs_fwd_deps( enums, solvmtx, rhsb, cblk, fcbk );
474 
475  /* Free the CBLK_RECV */
476  memFree_null( recvbuf );
477 }
478 
479 /**
480  *******************************************************************************
481  *
482  * @brief Handle a finished request.
483  *
484  * If cblktype & CBLK_FANIN : Will deallocate the coeftab
485  * If cblktype & CBLK_RECV : Will add cblk and deallocate the coeftab
486  *
487  *******************************************************************************
488  *
489  * @param[in] enums
490  * Enums needed for the solve.
491  *
492  * @param[inout] solvmtx
493  * The solver matrix structure.
494  *
495  * @param[in] rhsb
496  * The pointer to the rhs data structure that holds the vectors of the
497  * right hand side.
498  *
499  * @param[in] threadid
500  * Id of the thread calling this method.
501  *
502  * @param[in] outcount
503  * Amount of finshed requests.
504  *
505  * @param[in] indexes
506  * Array of completed requests.
507  *
508  * @param[in] statuses
509  * Array of statuses for the completed requests.
510  *
511  *******************************************************************************
512  *
513  * @retval The amount of finished requests updated.
514  *
515  *******************************************************************************/
516 static inline int
517 cpucblk_drequest_rhs_fwd_handle( const args_solve_t *enums,
518  SolverMatrix *solvmtx,
519  pastix_rhs_t rhsb,
520  int threadid,
521  int outcount,
522  const int *indexes,
523  const MPI_Status *statuses )
524 {
525  pastix_int_t i, reqid;
526  int nbrequest = outcount;
527 
528  for ( i = 0; i < outcount; i++ ) {
529  reqid = indexes[i];
530 
531  /*
532  * Handle the reception
533  */
534  if ( solvmtx->reqidx[reqid] == -1 ) {
535  /* We're on a cblk recv, copy datas and restart communications */
536  double *recvbuf;
537  MPI_Status status;
538  int size;
539 
540  memcpy( &status, statuses + i, sizeof(MPI_Status) );
541  MPI_Get_count( &status, PASTIX_MPI_DOUBLE, &size );
542 
543  MALLOC_INTERN( recvbuf, size, double );
544  memcpy( recvbuf, solvmtx->rcoeftab, size * sizeof(double) );
545 
546  solvmtx->recvcnt--;
547 
548  /* Let's restart the communication */
549  assert( solvmtx->recvcnt >= 0 );
550  if ( solvmtx->recvcnt > 0 ) {
551  MPI_Start( solvmtx->reqtab + reqid );
552  nbrequest--;
553  }
554  else {
555  MPI_Request_free( solvmtx->reqtab + reqid );
556  solvmtx->reqtab[reqid] = MPI_REQUEST_NULL;
557  }
558 
559  cpucblk_drequest_rhs_fwd_handle_recv( enums, solvmtx, rhsb,
560  threadid, &status, recvbuf );
561  }
562  /*
563  * Handle the emission
564  */
565  else {
566  SolverCblk *cblk = solvmtx->cblktab + solvmtx->reqidx[ reqid ];
567  assert( cblk->cblktype & CBLK_FANIN );
568 
569  cpucblk_drequest_rhs_fwd_handle_send( enums, solvmtx, rhsb, cblk );
570 
571 #if !defined(NDEBUG)
572  solvmtx->reqidx[ reqid ] = -1;
573 #endif
574  solvmtx->fanincnt--;
575  }
576  }
577 
578  return nbrequest;
579 }
580 
581 /**
582  *******************************************************************************
583  *
584  * @ingroup kernel_fact
585  * @brief Progress communications for one process
586  *
587  * If a communication is completed, it will be treated.
588  * If cblktype & CBLK_FANIN : Will deallocate coeftab
589  * If cblktype & CBLK_RECV : Will add cblk to fcblk
590  *
591  *******************************************************************************
592  *
593  * @param[in] enums
594  * Enums needed for the solve.
595  *
596  * @param[inout] solvmtx
597  * The solver matrix structure.
598  *
599  * @param[in] rhsb
600  * The pointer to the rhs data structure that holds the vectors of the
601  * right hand side.
602  *
603  * @param[in] threadid
604  * Id of the thread calling this method.
605  *
606  *******************************************************************************/
607 void
608 cpucblk_dmpi_rhs_fwd_progress( const args_solve_t *enums,
609  SolverMatrix *solvmtx,
610  pastix_rhs_t rhsb,
611  int threadid )
612 {
613  pthread_t tid = pthread_self();
614  int outcount = 1;
615  int nbrequest, nbfree;
616  int indexes[ solvmtx->reqnbr ];
617  MPI_Status statuses[ solvmtx->reqnbr ];
618 
619  /* Check if someone is already communicating or not */
620  pthread_mutex_lock( &pastix_comm_lock );
621  if ( pastix_comm_tid == (pthread_t)-1 ) {
622  pastix_comm_tid = tid;
623  }
624  pthread_mutex_unlock( &pastix_comm_lock );
625 
626  if ( tid != pastix_comm_tid ) {
627  pastix_yield();
628  return;
629  }
630 
631  /*
632  * Let's register the number of active requests.
633  * We now suppose that the current thread is working on the first nbrequest
634  * active in the reqtab array. Additional requests can be posted during this
635  * progression, but it will be with a larger index. Thus, we do not need to
636  * protect every changes in these requests.
637  * When this is done, the requests arrays is locked to be packed, and the
638  * number of requests is updated for the next round.
639  */
640  pastix_atomic_lock( &(solvmtx->reqlock) );
641  nbrequest = solvmtx->reqnum;
642  pastix_atomic_unlock( &(solvmtx->reqlock) );
643 
644  while( (outcount > 0) && (nbrequest > 0) )
645  {
646  MPI_Testsome( nbrequest, solvmtx->reqtab, &outcount, indexes, statuses );
647  nbfree = 0;
648 
649  /* Handle all the completed requests */
650  if ( outcount > 0 ) {
651  nbfree = cpucblk_drequest_rhs_fwd_handle( enums, solvmtx, rhsb, threadid,
652  outcount, indexes, statuses );
653  }
654 
655  /*
656  * Pack the request arrays, and update the number of active requests by
657  * removing the completed ones
658  */
659  pastix_atomic_lock( &(solvmtx->reqlock) );
660  if ( nbfree > 0 ) {
661  cpucblk_dupdate_reqtab( solvmtx );
662  }
663  nbrequest = solvmtx->reqnum;
664  pastix_atomic_unlock( &(solvmtx->reqlock) );
665  }
666 
667  pastix_comm_tid = (pthread_t)-1;
668  pastix_yield();
669 }
670 #endif /* defined(PASTIX_WITH_MPI) */
671 
672 /**
673  *******************************************************************************
674  *
675  * @brief Wait for incoming dependencies, and return when cblk->ctrbcnt has
676  * reached 0.
677  *
678  *******************************************************************************
679  *
680  * @param[in] rank
681  * The rank of the current thread.
682  *
683  * @param[in] enums
684  * Enums needed for the solve.
685  *
686  * @param[inout] solvmtx
687  * The solver matrix structure.
688  *
689  * @param[inout] cblk
690  * The column block that contribute to fcblk.
691  *
692  * @param[in] rhsb
693  * The pointer to the rhs data structure that holds the vectors of the
694  * right hand side.
695  *
696  *******************************************************************************
697  *
698  * @return 1 if the cblk is a fanin, 0 otherwise
699  *
700  *******************************************************************************/
701 int
703  const args_solve_t *enums,
704  SolverMatrix *solvmtx,
705  SolverCblk *cblk,
706  pastix_rhs_t rhsb )
707 {
708 #if defined(PASTIX_WITH_MPI)
709  if ( cblk->cblktype & CBLK_FANIN ) {
710  /*
711  * We are in the sequential case, we progress on communications and
712  * return if nothing.
713  */
714  //cpucblk_dtestsome( side, solvmtx );
715  return 1;
716  }
717 
718  if ( cblk->cblktype & CBLK_RECV ) {
719  return 1;
720  }
721 
722  /* Make sure we receive every contribution */
723  while( cblk->ctrbcnt > 0 ) {
724  cpucblk_dmpi_rhs_fwd_progress( enums, solvmtx, rhsb, rank );
725  }
726 #else
727  assert( !(cblk->cblktype & (CBLK_FANIN | CBLK_RECV)) );
728  do { pastix_yield(); } while( cblk->ctrbcnt > 0 );
729 #endif
730 
731  (void)rank;
732  (void)enums;
733  (void)solvmtx;
734  (void)rhsb;
735 
736  return 0;
737 }
738 
739 /**
740  *******************************************************************************
741  *
742  * @brief Release the dependencies of the given cblk after an update.
743  *
744  *******************************************************************************
745  *
746  * @param[in] enums
747  * Enums needed for the solve.
748  *
749  * @param[inout] solvmtx
750  * The solver matrix structure.
751  *
752  * @param[in] rhsb
753  * The pointer to the rhs data structure that holds the vectors of the
754  * right hand side.
755  *
756  * @param[in] cblk
757  * The column block that contribute to fcblk.
758  *
759  * @param[inout] fcbk
760  * The facing column block that is updated by cblk.
761  *
762  *******************************************************************************/
763 void
765  SolverMatrix *solvmtx,
766  pastix_rhs_t rhsb,
767  const SolverCblk *cblk,
768  SolverCblk *fcbk )
769 {
770  int32_t ctrbcnt;
771  ctrbcnt = pastix_atomic_dec_32b( &(fcbk->ctrbcnt) );
772  if ( !ctrbcnt ) {
773 #if defined(PASTIX_WITH_MPI)
774  if ( ( fcbk->cblktype & CBLK_FANIN ) &&
775  ( enums->solve_step == PastixSolveForward ) ) {
776  cpucblk_disend_rhs_fwd( solvmtx, rhsb, fcbk );
777  return;
778  }
779 #else
780  (void)enums;
781  (void)rhsb;
782 #endif
783  if ( solvmtx->computeQueue ) {
784  pastix_queue_t *queue = solvmtx->computeQueue[ cblk->threadid ];
785  assert( fcbk->priority != -1 );
786  pqueuePush1( queue, fcbk - solvmtx->cblktab, fcbk->priority );
787  }
788  }
789 }
790 
791 /**
792  *******************************************************************************
793  *
794  * @brief Waitall routine for current cblk request
795  *
796  * It may be possible that some cblk will not be deallocated with the static
797  * scheduler. So a cleanup may be necessary.
798  *
799  *******************************************************************************
800  *
801  * @param[in] enums
802  * Enums needed for the solve.
803  *
804  * @param[in] sched
805  * Define which sched is used
806  * @arg PastixSchedSequential if sequential
807  * @arg PastixSchedStatic if multi-threaded static scheduler
808  * @arg PastixSchedDynamic if multi-threaded dynamic scheduler
809  * No other scheduler is supported.
810  *
811  * @param[inout] solvmtx
812  * The solver matrix structure.
813  *
814  * @param[in] rhsb
815  * The pointer to the rhs data structure that holds the vectors of the
816  * right hand side.
817  *
818  *******************************************************************************/
819 void
821  pastix_int_t sched,
822  SolverMatrix *solvmtx,
823  pastix_rhs_t rhsb )
824 {
825  if ( (sched != PastixSchedSequential) &&
826  (sched != PastixSchedStatic) &&
827  (sched != PastixSchedDynamic) )
828  {
829  return;
830  }
831 #if defined(PASTIX_WITH_MPI)
832  pastix_int_t i;
833  int rc;
834  SolverCblk *cblk;
835  int reqnbr = solvmtx->reqnum;
836  MPI_Status status;
837 
838 #if defined(PASTIX_DEBUG_MPI)
839  fprintf( stderr, "[%2d] Wait for all pending communications\n",
840  solvmtx->clustnum );
841 #endif
842 
843  for ( i=0; i<reqnbr; i++ ) {
844  if ( solvmtx->reqtab[i] == MPI_REQUEST_NULL ) {
845  assert( 0 /* MPI_REQUEST_NULL should have been pushed to the end */ );
846  solvmtx->reqnum--;
847  continue;
848  }
849 
850  /* Make sure that we don't have an already cleaned request in dynamic */
851  assert( solvmtx->reqidx[i] != -1 );
852 
853  rc = MPI_Wait( solvmtx->reqtab + i, &status );
854  assert( rc == MPI_SUCCESS );
855 
856  cblk = solvmtx->cblktab + solvmtx->reqidx[i];
857 
858  /* We should wait only for fanin */
859  assert( cblk->cblktype & CBLK_FANIN );
860 
861  cpucblk_drequest_rhs_fwd_handle_send( enums, solvmtx, rhsb, cblk );
862 
863  solvmtx->reqnum--;
864  }
865  assert( solvmtx->reqnum == 0 );
866  (void)rc;
867 #else
868  (void)enums;
869  (void)solvmtx;
870  (void)rhsb;
871 #endif
872 }
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
static void pqueuePush1(pastix_queue_t *q, pastix_int_t elt, double key1)
Push an element with a single key.
Definition: queue.h:64
Queue structure.
Definition: queue.h:38
int core_dgeadd(pastix_trans_t trans, pastix_int_t M, pastix_int_t N, double alpha, const double *A, pastix_int_t LDA, double beta, double *B, pastix_int_t LDB)
Add two matrices together.
Definition: core_dgeadd.c:78
int cpucblk_dincoming_rhs_fwd_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 cpucblk_dsend_rhs_backward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Send the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_drecv_rhs_forward(const SolverMatrix *solvmtx, SolverCblk *cblk, double *work, pastix_rhs_t rhsb)
Receive the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_drecv_rhs_backward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Receive the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_dsend_rhs_forward(const SolverMatrix *solvmtx, SolverCblk *cblk, pastix_rhs_t rhsb)
Send the rhs associated to a cblk->lcolidx to the remote node.
void cpucblk_drelease_rhs_fwd_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_drequest_rhs_fwd_cleanup(const args_solve_t *enums, pastix_int_t sched, SolverMatrix *solvmtx, pastix_rhs_t rhsb)
Waitall routine for current cblk request.
@ PastixNoTrans
Definition: api.h:445
@ PastixSchedStatic
Definition: api.h:335
@ PastixSchedDynamic
Definition: api.h:338
@ PastixSchedSequential
Definition: api.h:334
void ** cblkb
Definition: pastixdata.h:157
pastix_int_t ld
Definition: pastixdata.h:155
pastix_int_t n
Definition: pastixdata.h:154
Main PaStiX RHS structure.
Definition: pastixdata.h:150
pastix_int_t reqnbr
Definition: solver.h:266
pastix_atomic_lock_t reqlock
Definition: solver.h:268
pastix_int_t priority
Definition: solver.h:177
static pastix_int_t cblk_colnbr(const SolverCblk *cblk)
Compute the number of columns in a column block.
Definition: solver.h:324
pastix_int_t fcblknm
Definition: solver.h:140
MPI_Request * reqtab
Definition: solver.h:264
pastix_int_t recvcnt
Definition: solver.h:213
void * rcoeftab
Definition: solver.h:269
pastix_int_t * gcbl2loc
Definition: solver.h:228
pastix_int_t fanincnt
Definition: solver.h:210
volatile int32_t ctrbcnt
Definition: solver.h:158
SolverBlok * fblokptr
Definition: solver.h:163
volatile int32_t reqnum
Definition: solver.h:267
pastix_int_t gcblknum
Definition: solver.h:169
int threadid
Definition: solver.h:176
pastix_int_t lcolidx
Definition: solver.h:165
pastix_int_t bcscnum
Definition: solver.h:170
SolverCblk *restrict cblktab
Definition: solver.h:222
int8_t cblktype
Definition: solver.h:159
pastix_int_t * reqidx
Definition: solver.h:265
int ownerid
Definition: solver.h:175
Arguments for the solve.
Definition: solver.h:85
Solver column block structure.
Definition: solver.h:156
Solver column block structure.
Definition: solver.h:200