xref: /petsc/src/ksp/pc/impls/tfs/comm.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
1 
2 /***********************************comm.c*************************************
3 
4 Author: Henry M. Tufo III
5 
6 e-mail: hmt@cs.brown.edu
7 
8 snail-mail:
9 Division of Applied Mathematics
10 Brown University
11 Providence, RI 02912
12 
13 Last Modification:
14 11.21.97
15 ***********************************comm.c*************************************/
16 #include <../src/ksp/pc/impls/tfs/tfs.h>
17 
18 
19 /* global program control variables - explicitly exported */
20 PetscMPIInt PCTFS_my_id            = 0;
21 PetscMPIInt PCTFS_num_nodes        = 1;
22 PetscMPIInt PCTFS_floor_num_nodes  = 0;
23 PetscMPIInt PCTFS_i_log2_num_nodes = 0;
24 
25 /* global program control variables */
26 static PetscInt p_init = 0;
27 static PetscInt modfl_num_nodes;
28 static PetscInt edge_not_pow_2;
29 
30 static PetscInt edge_node[sizeof(PetscInt)*32];
31 
32 /***********************************comm.c*************************************/
33 PetscErrorCode PCTFS_comm_init(void)
34 {
35 
36   if (p_init++) PetscFunctionReturn(0);
37 
38   MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes);
39   MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id);
40 
41   if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");
42 
43   PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);
44 
45   PCTFS_floor_num_nodes  = 1;
46   PCTFS_i_log2_num_nodes = modfl_num_nodes = 0;
47   while (PCTFS_floor_num_nodes <= PCTFS_num_nodes) {
48     edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes;
49     PCTFS_floor_num_nodes           <<= 1;
50     PCTFS_i_log2_num_nodes++;
51   }
52 
53   PCTFS_i_log2_num_nodes--;
54   PCTFS_floor_num_nodes >>= 1;
55   modfl_num_nodes         = (PCTFS_num_nodes - PCTFS_floor_num_nodes);
56 
57   if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes)) edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1);
58   else if (PCTFS_my_id >= PCTFS_floor_num_nodes) edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1);
59   else edge_not_pow_2 = 0;
60   PetscFunctionReturn(0);
61 }
62 
63 /***********************************comm.c*************************************/
64 PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
65 {
66   PetscInt   mask, edge;
67   PetscInt   type, dest;
68   vfp        fp;
69   MPI_Status status;
70   PetscInt   ierr;
71 
72   PetscFunctionBegin;
73   /* ok ... should have some data, work, and operator(s) */
74   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
75 
76   /* non-uniform should have at least two entries */
77   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?");
78 
79   /* check to make sure comm package has been initialized */
80   if (!p_init) PCTFS_comm_init();
81 
82   /* if there's nothing to do return */
83   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
84 
85   /* a negative number if items to send ==> fatal */
86   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n);
87 
88   /* advance to list of n operations for custom */
89   if ((type=oprs[0])==NON_UNIFORM) oprs++;
90 
91   /* major league hack */
92   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: Could not retrieve function pointer!\n");
93 
94   /* all msgs will be of the same length */
95   /* if not a hypercube must colapse partial dim */
96   if (edge_not_pow_2) {
97     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
98       ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
99     } else {
100       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
101       (*fp)(vals,work,n,oprs);
102     }
103   }
104 
105   /* implement the mesh fan in/out exchange algorithm */
106   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
107     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
108       dest = PCTFS_my_id^mask;
109       if (PCTFS_my_id > dest) {
110         ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
111       } else {
112         ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
113         (*fp)(vals, work, n, oprs);
114       }
115     }
116 
117     mask=PCTFS_floor_num_nodes>>1;
118     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
119       if (PCTFS_my_id%mask) continue;
120 
121       dest = PCTFS_my_id^mask;
122       if (PCTFS_my_id < dest) {
123         ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
124       } else {
125         ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
126       }
127     }
128   }
129 
130   /* if not a hypercube must expand to partial dim */
131   if (edge_not_pow_2) {
132     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
133       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
134     } else {
135       ierr = MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
136     }
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 /***********************************comm.c*************************************/
142 PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
143 {
144   PetscInt       mask, edge;
145   PetscInt       type, dest;
146   vfp            fp;
147   MPI_Status     status;
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   /* ok ... should have some data, work, and operator(s) */
152   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
153 
154   /* non-uniform should have at least two entries */
155   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?");
156 
157   /* check to make sure comm package has been initialized */
158   if (!p_init) PCTFS_comm_init();
159 
160   /* if there's nothing to do return */
161   if ((PCTFS_num_nodes<2)||(!n)) PetscFunctionReturn(0);
162 
163   /* a negative number of items to send ==> fatal */
164   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);
165 
166   /* advance to list of n operations for custom */
167   if ((type=oprs[0])==NON_UNIFORM) oprs++;
168 
169   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: Could not retrieve function pointer!\n");
170 
171   /* all msgs will be of the same length */
172   /* if not a hypercube must colapse partial dim */
173   if (edge_not_pow_2) {
174     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
175       ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
176     } else {
177       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
178       (*fp)(vals,work,n,oprs);
179     }
180   }
181 
182   /* implement the mesh fan in/out exchange algorithm */
183   if (PCTFS_my_id<PCTFS_floor_num_nodes) {
184     for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1) {
185       dest = PCTFS_my_id^mask;
186       if (PCTFS_my_id > dest) {
187         ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
188       } else {
189         ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
190         (*fp)(vals, work, n, oprs);
191       }
192     }
193 
194     mask=PCTFS_floor_num_nodes>>1;
195     for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1) {
196       if (PCTFS_my_id%mask) continue;
197 
198       dest = PCTFS_my_id^mask;
199       if (PCTFS_my_id < dest) {
200         ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
201       } else {
202         ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
203       }
204     }
205   }
206 
207   /* if not a hypercube must expand to partial dim */
208   if (edge_not_pow_2) {
209     if (PCTFS_my_id >= PCTFS_floor_num_nodes) {
210       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);CHKERRQ(ierr);
211     } else {
212       ierr = MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
213     }
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 /***********************************comm.c*************************************/
219 PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
220 {
221   PetscInt       mask, edge;
222   PetscInt       type, dest;
223   vfp            fp;
224   MPI_Status     status;
225   PetscErrorCode ierr;
226 
227   PetscFunctionBegin;
228   /* ok ... should have some data, work, and operator(s) */
229   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
230 
231   /* non-uniform should have at least two entries */
232   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?");
233 
234   /* check to make sure comm package has been initialized */
235   if (!p_init) PCTFS_comm_init();
236 
237   /* if there's nothing to do return */
238   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
239 
240   /* the error msg says it all!!! */
241   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");
242 
243   /* a negative number of items to send ==> fatal */
244   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n);
245 
246   /* can't do more dimensions then exist */
247   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
248 
249   /* advance to list of n operations for custom */
250   if ((type=oprs[0])==NON_UNIFORM) oprs++;
251 
252   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: Could not retrieve function pointer!\n");
253 
254   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
255     dest = PCTFS_my_id^mask;
256     if (PCTFS_my_id > dest) {
257       ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
258     } else {
259       ierr = MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
260       (*fp)(vals, work, n, oprs);
261     }
262   }
263 
264   if (edge==dim) mask>>=1;
265   else {
266     while (++edge<dim) mask<<=1;
267   }
268 
269   for (edge=0; edge<dim; edge++,mask>>=1) {
270     if (PCTFS_my_id%mask) continue;
271 
272     dest = PCTFS_my_id^mask;
273     if (PCTFS_my_id < dest) {
274       ierr = MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
275     } else {
276       ierr = MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
277     }
278   }
279   PetscFunctionReturn(0);
280 }
281 
282 /******************************************************************************/
283 PetscErrorCode PCTFS_ssgl_radd(PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
284 {
285   PetscInt       edge, type, dest, mask;
286   PetscInt       stage_n;
287   MPI_Status     status;
288   PetscErrorCode ierr;
289 
290   PetscFunctionBegin;
291   /* check to make sure comm package has been initialized */
292   if (!p_init) PCTFS_comm_init();
293 
294   /* all msgs are *NOT* the same length */
295   /* implement the mesh fan in/out exchange algorithm */
296   for (mask=0, edge=0; edge<level; edge++, mask++) {
297     stage_n = (segs[level] - segs[edge]);
298     if (stage_n && !(PCTFS_my_id & mask)) {
299       dest = edge_node[edge];
300       type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge);
301       if (PCTFS_my_id>dest) {
302         ierr = MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);CHKERRQ(ierr);
303       } else {
304         type =  type - PCTFS_my_id + dest;
305         ierr = MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
306         PCTFS_rvec_add(vals+segs[edge], work, stage_n);
307       }
308     }
309     mask <<= 1;
310   }
311   mask>>=1;
312   for (edge=0; edge<level; edge++) {
313     stage_n = (segs[level] - segs[level-1-edge]);
314     if (stage_n && !(PCTFS_my_id & mask)) {
315       dest = edge_node[level-edge-1];
316       type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge);
317       if (PCTFS_my_id<dest) {
318         ierr = MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);CHKERRQ(ierr);
319       } else {
320         type =  type - PCTFS_my_id + dest;
321         ierr = MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
322       }
323     }
324     mask >>= 1;
325   }
326   PetscFunctionReturn(0);
327 }
328 
329 /***********************************comm.c*************************************/
330 PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
331 {
332   PetscInt       mask, edge;
333   PetscInt       type, dest;
334   vfp            fp;
335   MPI_Status     status;
336   PetscErrorCode ierr;
337 
338   PetscFunctionBegin;
339   /* ok ... should have some data, work, and operator(s) */
340   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);
341 
342   /* non-uniform should have at least two entries */
343   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?");
344 
345   /* check to make sure comm package has been initialized */
346   if (!p_init) PCTFS_comm_init();
347 
348   /* if there's nothing to do return */
349   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0)) PetscFunctionReturn(0);
350 
351   /* the error msg says it all!!! */
352   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");
353 
354   /* a negative number of items to send ==> fatal */
355   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n);
356 
357   /* can't do more dimensions then exist */
358   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);
359 
360   /* advance to list of n operations for custom */
361   if ((type=oprs[0])==NON_UNIFORM) oprs++;
362 
363   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: Could not retrieve function pointer!\n");
364 
365   for (mask=1,edge=0; edge<dim; edge++,mask<<=1) {
366     dest = PCTFS_my_id^mask;
367     if (PCTFS_my_id > dest) {
368       ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
369     } else {
370       ierr = MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);CHKERRQ(ierr);
371       (*fp)(vals, work, n, oprs);
372     }
373   }
374 
375   if (edge==dim) mask>>=1;
376   else {
377     while (++edge<dim) mask<<=1;
378   }
379 
380   for (edge=0; edge<dim; edge++,mask>>=1) {
381     if (PCTFS_my_id%mask) continue;
382 
383     dest = PCTFS_my_id^mask;
384     if (PCTFS_my_id < dest) {
385       ierr = MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);CHKERRQ(ierr);
386     } else {
387       ierr = MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);CHKERRQ(ierr);
388     }
389   }
390   PetscFunctionReturn(0);
391 }
392