PaStiX Handbook  6.3.2
solver_recv.c
Go to the documentation of this file.
1 /**
2  * @file solver_recv.c
3  *
4  * PaStiX solver reception structure management.
5  *
6  * @copyright 2004-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
7  * Univ. Bordeaux. All rights reserved.
8  *
9  * @version 6.3.2
10  * @author Mathieu Faverge
11  * @author Tony Delarue
12  * @author Alycia Lisito
13  * @date 2023-12-01
14  *
15  * @addtogroup blend_dev_solver
16  * @{
17  *
18  **/
19 #include "common.h"
20 #include "symbol/symbol.h"
22 
23 /**
24  *******************************************************************************
25  *
26  * @brief Update columns indices of a reception/fanin cblk.
27  *
28  *******************************************************************************
29  *
30  * @param[inout] cblk
31  * The fanin or recv cblk to update. Must have been initialized first.
32  *
33  * @param[in] fcolnum
34  * The first column index of the contribution block.
35  *
36  * @param[in] lcolnum
37  * The last column index of the contribution block.
38  *
39  *******************************************************************************/
40 static inline void
42  pastix_int_t fcolnum,
43  pastix_int_t lcolnum )
44 {
45  cblk->fcolnum = pastix_imin( fcolnum, cblk->fcolnum );
46  cblk->lcolnum = pastix_imax( lcolnum, cblk->lcolnum );
47 }
48 
49 /**
50  *******************************************************************************
51  *
52  * @brief Update rows indices of a reception/fanin blok.
53  *
54  *******************************************************************************
55  *
56  * @param[inout] blok
57  * The fanin or recv blok to update. Must have been initialized first.
58  *
59  * @param[in] frownum
60  * The first row index of the contribution block.
61  *
62  * @param[in] lrownum
63  * The last row index of the contribution block.
64  *
65  *******************************************************************************/
66 static inline void
68  pastix_int_t frownum,
69  pastix_int_t lrownum )
70 {
71  blok->frownum = pastix_imin( frownum, blok->frownum );
72  blok->lrownum = pastix_imax( lrownum, blok->lrownum );
73 }
74 
75 /**
76  *******************************************************************************
77  *
78  * @brief Create a new reception/fanin cblk and initialize to the default values.
79  *
80  *******************************************************************************
81  *
82  * @param[in] symbmtx
83  * The symbol matrix pointer (used to access the bloktab array)
84  *
85  * @param[in] cblk
86  * The symbol cblk used as a template to create the fanin/recv cblk.
87  *
88  *******************************************************************************
89  *
90  * @return The pointer to the recv/fanin cblk initialized with respect to the
91  * input symbcblk.
92  *
93  *******************************************************************************/
94 static inline solver_cblk_recv_t *
96  const symbol_cblk_t *cblk )
97 {
98  solver_cblk_recv_t *cblkrecv;
99  symbol_blok_t *blok;
100  pastix_int_t bloknbr = cblk[1].bloknum - cblk[0].bloknum;
101  pastix_int_t i;
102 
103  assert( bloknbr >= 1 );
104  cblkrecv = malloc( sizeof(solver_cblk_recv_t) + ((bloknbr-1) * sizeof(solver_blok_recv_t)) );
105  cblkrecv->next = NULL;
106  cblkrecv->ownerid = -1;
107  cblkrecv->fcolnum = cblk->lcolnum + 1;
108  cblkrecv->lcolnum = cblk->fcolnum - 1;
109 
110  blok = symbmtx->bloktab + cblk->bloknum;
111  for( i=0; i<bloknbr; i++, blok++ ) {
112  cblkrecv->bloktab[i].frownum = blok->lrownum + 1;
113  cblkrecv->bloktab[i].lrownum = blok->frownum - 1;
114  }
115 
116  return cblkrecv;
117 }
118 
119 /**
120  *******************************************************************************
121  *
122  * @brief TODO
123  *
124  *******************************************************************************
125  *
126  * @param[in] rcblk
127  * TODO
128  *
129  * @param[in] symbmtx
130  * TODO
131  *
132  * @param[in] cblk
133  * TODO
134  *
135  * @param[in] blok
136  * TODO
137  *
138  * @param[in] fcblk
139  * TODO
140  *
141  *******************************************************************************/
142 static inline void
144  const symbol_matrix_t *symbmtx,
145  const symbol_cblk_t *cblk,
146  const symbol_blok_t *blok,
147  const symbol_cblk_t *fcblk )
148 {
149  symbol_blok_t *lblok, *fblok;
150  pastix_int_t i;
151 
152  /* Let's update the columns of the cblk */
153  solver_recv_update_cols( rcblk, blok->frownum, blok->lrownum );
154 
155  /* Let's iterate on the blocks to update the rows */
156  lblok = symbmtx->bloktab + cblk[1].bloknum;
157  fblok = symbmtx->bloktab + fcblk->bloknum;
158  i = 0;
159  for( ; blok<lblok; blok++ ) {
160  /* Look for the facing blok in the symbol */
161  while(!is_symbblock_inside_fblock( blok, fblok )) {
162  i++;
163  fblok++;
164  assert( fcblk[0].bloknum + i < fcblk[1].bloknum );
165  }
166 
167  /* Update the rows in the recv block */
168  solver_recv_update_rows( &(rcblk->bloktab[i]), blok->frownum, blok->lrownum );
169  }
170 }
171 
172 /**
173  *******************************************************************************
174  *
175  * @brief Register a new contribution to a fanin cblk
176  *
177  *******************************************************************************
178  *
179  * @param[inout] faninptr
180  * The location of the pointer of the fanin cblk.
181  * On entry, points to NULL if the fanin does not exist, to the
182  * existing fanin otherwise.
183  * On exit, points to the updated fanin cblk.
184  *
185  * @param[in] symbmtx
186  * The symbol matrix pointer (used to access the bloktab array)
187  *
188  * @param[in] cblk
189  * The symbol cblk that holds the block responsible for the contribution.
190  *
191  * @param[in] blok
192  * The symbol blok responsible for the contribution.
193  *
194  * @param[in] fcblk
195  * The facing symbol cblk that corresponds to the faninptr, and that is
196  * updated by the contribution.
197  *
198  * @param[in] ownerid
199  * The index of the cluster that owns the fcblk.
200  *
201  *******************************************************************************/
202 void
204  const symbol_matrix_t *symbmtx,
205  const symbol_cblk_t *cblk,
206  const symbol_blok_t *blok,
207  const symbol_cblk_t *fcblk,
208  int ownerid )
209 {
210  if ( *faninptr == NULL ) {
211  *faninptr = solver_recv_cblk_init( symbmtx, fcblk );
212  (*faninptr)->ownerid = ownerid;
213  }
214  assert( (*faninptr)->ownerid == ownerid );
215 
216  solver_recv_add_contrib( *faninptr, symbmtx, cblk, blok, fcblk );
217 }
218 
219 /**
220  *******************************************************************************
221  *
222  * @brief Register a new contribution to a recv cblk
223  *
224  *******************************************************************************
225  *
226  * @param[inout] recvptr
227  * The location of the pointer of the recv cblks list.
228  * On entry, points to NULL if no recv cblk has been registered yet, to the
229  * head of the list otherwise.
230  * On exit, points to the updated list of recv cblks.
231  *
232  * @param[in] symbmtx
233  * The symbol matrix pointer (used to access the bloktab array)
234  *
235  * @param[in] cblk
236  * The symbol cblk that holds the block responsible for the contribution.
237  *
238  * @param[in] blok
239  * The symbol blok responsible for the contribution.
240  *
241  * @param[in] fcblk
242  * The facing symbol cblk that corresponds to the recvptr, and that is
243  * updated by the contribution.
244  *
245  * @param[in] ownerid
246  * The index of the cluster that owns the original contribution (cblk).
247  *
248  *******************************************************************************/
249 void
251  const symbol_matrix_t *symbmtx,
252  const symbol_cblk_t *cblk,
253  const symbol_blok_t *blok,
254  const symbol_cblk_t *fcblk,
255  int ownerid )
256 {
257  solver_cblk_recv_t *prev, *rcblk;
258 
259  prev = NULL;
260  rcblk = *recvptr;
261  while( (rcblk != NULL) && (rcblk->ownerid != ownerid) ) {
262  prev = rcblk;
263  rcblk = rcblk->next;
264  }
265 
266  /* Create a new cblk if not found */
267  if ( rcblk == NULL ) {
268  rcblk = solver_recv_cblk_init( symbmtx, fcblk );
269  rcblk->ownerid = ownerid;
270  if ( prev == NULL ) {
271  /* Head of the list */
272  *recvptr = rcblk;
273  }
274  else {
275  assert( prev->next == NULL );
276  prev->next = rcblk;
277  }
278  }
279 
280  assert( rcblk->ownerid == ownerid );
281  solver_recv_add_contrib( rcblk, symbmtx, cblk, blok, fcblk );
282 }
283 
284 /**
285  *******************************************************************************
286  *
287  * @brief Compute the number of valid blocks in fanin/recv cblk
288  *
289  *******************************************************************************
290  *
291  * @param[in] ftgtptr
292  * The pointer to the fan-in contribution to study.
293  *
294  * @param[in] symbcblk
295  * The original symbol cblk associated to the fan-in
296  *
297  * @param[in] symbblok
298  * The first symbol blok of the symbcblk.
299  *
300  *******************************************************************************
301  *
302  * @return The number of non empty blocks in the fan-in.
303  *
304  *******************************************************************************/
305 int
307  const symbol_cblk_t *symbcblk,
308  const symbol_blok_t *symbblok )
309 {
310  const solver_blok_recv_t *ftgtblok = ftgtptr->bloktab;
311  pastix_int_t j;
312  int bloknbr = 0;
313 
314  for ( j=symbcblk[0].bloknum; j<symbcblk[1].bloknum;
315  j++, symbblok++, ftgtblok++ )
316  {
317  if ( (ftgtblok->frownum <= ftgtblok->lrownum) &&
318  (ftgtblok->frownum >= symbblok->frownum) &&
319  (ftgtblok->lrownum <= symbblok->lrownum) )
320  {
321  bloknbr++;
322  }
323  }
324  assert( bloknbr >= 1 );
325  return bloknbr;
326 }
327 
328 /**
329  * @}
330  */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
int solver_recv_get_bloknbr(const solver_cblk_recv_t *ftgtptr, const symbol_cblk_t *symbcblk, const symbol_blok_t *symbblok)
Compute the number of valid blocks in fanin/recv cblk.
Definition: solver_recv.c:306
static void solver_recv_add_contrib(solver_cblk_recv_t *rcblk, const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk, const symbol_blok_t *blok, const symbol_cblk_t *fcblk)
TODO.
Definition: solver_recv.c:143
void solver_recv_update_fanin(solver_cblk_recv_t **faninptr, const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk, const symbol_blok_t *blok, const symbol_cblk_t *fcblk, int ownerid)
Register a new contribution to a fanin cblk.
Definition: solver_recv.c:203
void solver_recv_update_recv(solver_cblk_recv_t **recvptr, const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk, const symbol_blok_t *blok, const symbol_cblk_t *fcblk, int ownerid)
Register a new contribution to a recv cblk.
Definition: solver_recv.c:250
static void solver_recv_update_rows(solver_blok_recv_t *blok, pastix_int_t frownum, pastix_int_t lrownum)
Update rows indices of a reception/fanin blok.
Definition: solver_recv.c:67
static void solver_recv_update_cols(solver_cblk_recv_t *cblk, pastix_int_t fcolnum, pastix_int_t lcolnum)
Update columns indices of a reception/fanin cblk.
Definition: solver_recv.c:41
static solver_cblk_recv_t * solver_recv_cblk_init(const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk)
Create a new reception/fanin cblk and initialize to the default values.
Definition: solver_recv.c:95
pastix_int_t frownum
Definition: symbol.h:60
pastix_int_t lrownum
Definition: symbol.h:61
pastix_int_t bloknum
Definition: symbol.h:48
pastix_int_t fcolnum
Definition: symbol.h:46
pastix_int_t lcolnum
Definition: symbol.h:47
symbol_blok_t * bloktab
Definition: symbol.h:84
static int is_symbblock_inside_fblock(const symbol_blok_t *blok, const symbol_blok_t *fblok)
Check if a block is included inside another one.
Definition: symbol.h:185
Symbol block structure.
Definition: symbol.h:59
Symbol column block structure.
Definition: symbol.h:45
Symbol matrix structure.
Definition: symbol.h:77
pastix_int_t lrownum
Definition: solver.h:99
pastix_int_t frownum
Definition: solver.h:98
pastix_int_t fcolnum
Definition: solver.h:110
pastix_int_t lcolnum
Definition: solver.h:111
solver_blok_recv_t bloktab[1]
Definition: solver.h:112
Solver recv block structure.
Definition: solver.h:97
Solver recv column block structure.
Definition: solver.h:107