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