Actual source code: comm.c

petsc-3.3-p7 2013-05-11
  2: /***********************************comm.c*************************************

  4: Author: Henry M. Tufo III

  6: e-mail: hmt@cs.brown.edu

  8: snail-mail:
  9: Division of Applied Mathematics
 10: Brown University
 11: Providence, RI 02912

 13: Last Modification: 
 14: 11.21.97
 15: ***********************************comm.c*************************************/
 16: #include <../src/ksp/pc/impls/tfs/tfs.h>


 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;

 25: /* global program control variables */
 26: static PetscInt p_init = 0;
 27: static PetscInt modfl_num_nodes;
 28: static PetscInt edge_not_pow_2;

 30: static PetscInt edge_node[sizeof(PetscInt)*32];

 32: /***********************************comm.c*************************************/
 33: PetscErrorCode PCTFS_comm_init (void)
 34: {

 36:   if (p_init++)   return(0);

 38:   MPI_Comm_size(MPI_COMM_WORLD,&PCTFS_num_nodes);
 39:   MPI_Comm_rank(MPI_COMM_WORLD,&PCTFS_my_id);

 41:   if (PCTFS_num_nodes> (INT_MAX >> 1)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Can't have more then MAX_INT/2 nodes!!!");

 43:   PCTFS_ivec_zero((PetscInt*)edge_node,sizeof(PetscInt)*32);

 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:     {
 49:       edge_node[PCTFS_i_log2_num_nodes] = PCTFS_my_id ^ PCTFS_floor_num_nodes;
 50:       PCTFS_floor_num_nodes <<= 1;
 51:       PCTFS_i_log2_num_nodes++;
 52:     }

 54:   PCTFS_i_log2_num_nodes--;
 55:   PCTFS_floor_num_nodes >>= 1;
 56:   modfl_num_nodes = (PCTFS_num_nodes - PCTFS_floor_num_nodes);

 58:   if ((PCTFS_my_id > 0) && (PCTFS_my_id <= modfl_num_nodes))
 59:     {edge_not_pow_2=((PCTFS_my_id|PCTFS_floor_num_nodes)-1);}
 60:   else if (PCTFS_my_id >= PCTFS_floor_num_nodes)
 61:     {edge_not_pow_2=((PCTFS_my_id^PCTFS_floor_num_nodes)+1);
 62:     }
 63:   else
 64:     {edge_not_pow_2 = 0;}
 65:   return(0);
 66: }

 68: /***********************************comm.c*************************************/
 69: PetscErrorCode PCTFS_giop(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs)
 70: {
 71:   PetscInt   mask, edge;
 72:   PetscInt    type, dest;
 73:   vfp         fp;
 74:   MPI_Status  status;
 75:   PetscInt    ierr;

 78:   /* ok ... should have some data, work, and operator(s) */
 79:   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);

 81:   /* non-uniform should have at least two entries */
 82:   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: non_uniform and n=0,1?");

 84:   /* check to make sure comm package has been initialized */
 85:   if (!p_init)
 86:     {PCTFS_comm_init();}

 88:   /* if there's nothing to do return */
 89:   if ((PCTFS_num_nodes<2)||(!n))
 90:     {
 91:         return(0);
 92:     }


 95:   /* a negative number if items to send ==> fatal */
 96:   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop() :: n=%D<0?",n);

 98:   /* advance to list of n operations for custom */
 99:   if ((type=oprs[0])==NON_UNIFORM)
100:     {oprs++;}

102:   /* major league hack */
103:   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))) {
104:     PetscInfo(0,"PCTFS_giop() :: hope you passed in a rbfp!\n");
105:     fp = (vfp) oprs;
106:   }

108:   /* all msgs will be of the same length */
109:   /* if not a hypercube must colapse partial dim */
110:   if (edge_not_pow_2)
111:     {
112:       if (PCTFS_my_id >= PCTFS_floor_num_nodes)
113:         {MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);}
114:       else
115:         {
116:           MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2, MPI_COMM_WORLD,&status);
117:           (*fp)(vals,work,n,oprs);
118:         }
119:     }

121:   /* implement the mesh fan in/out exchange algorithm */
122:   if (PCTFS_my_id<PCTFS_floor_num_nodes)
123:     {
124:       for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1)
125:         {
126:           dest = PCTFS_my_id^mask;
127:           if (PCTFS_my_id > dest)
128:             {MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);}
129:           else
130:             {
131:               MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);
132:               (*fp)(vals, work, n, oprs);
133:             }
134:         }

136:       mask=PCTFS_floor_num_nodes>>1;
137:       for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1)
138:         {
139:           if (PCTFS_my_id%mask)
140:             {continue;}
141: 
142:           dest = PCTFS_my_id^mask;
143:           if (PCTFS_my_id < dest)
144:             {MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);}
145:           else
146:             {
147:               MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);
148:             }
149:         }
150:     }

152:   /* if not a hypercube must expand to partial dim */
153:   if (edge_not_pow_2)
154:     {
155:       if (PCTFS_my_id >= PCTFS_floor_num_nodes)
156:         {
157:           MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2,MPI_COMM_WORLD,&status);
158:         }
159:       else
160:         {MPI_Send(vals,n,MPIU_INT,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);}
161:     }
162:         return(0);
163: }

165: /***********************************comm.c*************************************/
166: PetscErrorCode PCTFS_grop(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs)
167: {
168:   PetscInt       mask, edge;
169:   PetscInt       type, dest;
170:   vfp            fp;
171:   MPI_Status     status;

175:   /* ok ... should have some data, work, and operator(s) */
176:   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);

178:   /* non-uniform should have at least two entries */
179:   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop() :: non_uniform and n=0,1?");

181:   /* check to make sure comm package has been initialized */
182:   if (!p_init)
183:     {PCTFS_comm_init();}

185:   /* if there's nothing to do return */
186:   if ((PCTFS_num_nodes<2)||(!n))
187:     {        return(0);}

189:   /* a negative number of items to send ==> fatal */
190:   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gdop() :: n=%D<0?",n);

192:   /* advance to list of n operations for custom */
193:   if ((type=oprs[0])==NON_UNIFORM)
194:     {oprs++;}

196:   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) {
197:     PetscInfo(0,"PCTFS_grop() :: hope you passed in a rbfp!\n");
198:     fp = (vfp) oprs;
199:   }

201:   /* all msgs will be of the same length */
202:   /* if not a hypercube must colapse partial dim */
203:   if (edge_not_pow_2)
204:     {
205:       if (PCTFS_my_id >= PCTFS_floor_num_nodes)
206:         {MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG0+PCTFS_my_id,MPI_COMM_WORLD);}
207:       else
208:         {
209:           MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG0+edge_not_pow_2,MPI_COMM_WORLD,&status);
210:           (*fp)(vals,work,n,oprs);
211:         }
212:     }

214:   /* implement the mesh fan in/out exchange algorithm */
215:   if (PCTFS_my_id<PCTFS_floor_num_nodes)
216:     {
217:       for (mask=1,edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask<<=1)
218:         {
219:           dest = PCTFS_my_id^mask;
220:           if (PCTFS_my_id > dest)
221:             {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);}
222:           else
223:             {
224:               MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);
225:               (*fp)(vals, work, n, oprs);
226:             }
227:         }

229:       mask=PCTFS_floor_num_nodes>>1;
230:       for (edge=0; edge<PCTFS_i_log2_num_nodes; edge++,mask>>=1)
231:         {
232:           if (PCTFS_my_id%mask)
233:             {continue;}
234: 
235:           dest = PCTFS_my_id^mask;
236:           if (PCTFS_my_id < dest)
237:             {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);}
238:           else
239:             {
240:               MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD, &status);
241:             }
242:         }
243:     }

245:   /* if not a hypercube must expand to partial dim */
246:   if (edge_not_pow_2)
247:     {
248:       if (PCTFS_my_id >= PCTFS_floor_num_nodes)
249:         {
250:           MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG5+edge_not_pow_2, MPI_COMM_WORLD,&status);
251:         }
252:       else
253:         {MPI_Send(vals,n,MPIU_SCALAR,edge_not_pow_2,MSGTAG5+PCTFS_my_id,MPI_COMM_WORLD);}
254:     }
255:         return(0);
256: }

258: /***********************************comm.c*************************************/
259: PetscErrorCode PCTFS_grop_hc(PetscScalar *vals, PetscScalar *work, PetscInt n, PetscInt *oprs, PetscInt dim)
260: {
261:   PetscInt       mask, edge;
262:   PetscInt       type, dest;
263:   vfp            fp;
264:   MPI_Status     status;

268:   /* ok ... should have some data, work, and operator(s) */
269:   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);

271:   /* non-uniform should have at least two entries */
272:   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: non_uniform and n=0,1?");

274:   /* check to make sure comm package has been initialized */
275:   if (!p_init)
276:     {PCTFS_comm_init();}

278:   /* if there's nothing to do return */
279:   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0))
280:     {return(0);}

282:   /* the error msg says it all!!! */
283:   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: PCTFS_num_nodes not a power of 2!?!");

285:   /* a negative number of items to send ==> fatal */
286:   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_grop_hc() :: n=%D<0?",n);

288:   /* can't do more dimensions then exist */
289:   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);

291:   /* advance to list of n operations for custom */
292:   if ((type=oprs[0])==NON_UNIFORM)
293:     {oprs++;}

295:   if (!(fp = (vfp) PCTFS_rvec_fct_addr(type))) {
296:     PetscInfo(0,"PCTFS_grop_hc() :: hope you passed in a rbfp!\n");
297:     fp = (vfp) oprs;
298:   }

300:   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
301:     {
302:       dest = PCTFS_my_id^mask;
303:       if (PCTFS_my_id > dest)
304:         {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);}
305:       else
306:         {
307:           MPI_Recv(work,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD,&status);
308:           (*fp)(vals, work, n, oprs);
309:         }
310:     }

312:   if (edge==dim)
313:     {mask>>=1;}
314:   else
315:     {while (++edge<dim) {mask<<=1;}}

317:   for (edge=0; edge<dim; edge++,mask>>=1)
318:     {
319:       if (PCTFS_my_id%mask)
320:         {continue;}
321: 
322:       dest = PCTFS_my_id^mask;
323:       if (PCTFS_my_id < dest)
324:         {MPI_Send(vals,n,MPIU_SCALAR,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);}
325:       else
326:         {
327:           MPI_Recv(vals,n,MPIU_SCALAR,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);
328:         }
329:     }
330:         return(0);
331: }

333: /******************************************************************************/
334: PetscErrorCode PCTFS_ssgl_radd( PetscScalar *vals,  PetscScalar *work,  PetscInt level, PetscInt *segs)
335: {
336:   PetscInt       edge, type, dest, mask;
337:   PetscInt       stage_n;
338:   MPI_Status     status;

342:   /* check to make sure comm package has been initialized */
343:   if (!p_init)
344:     {PCTFS_comm_init();}


347:   /* all msgs are *NOT* the same length */
348:   /* implement the mesh fan in/out exchange algorithm */
349:   for (mask=0, edge=0; edge<level; edge++, mask++)
350:     {
351:       stage_n = (segs[level] - segs[edge]);
352:       if (stage_n && !(PCTFS_my_id & mask))
353:         {
354:           dest = edge_node[edge];
355:           type = MSGTAG3 + PCTFS_my_id + (PCTFS_num_nodes*edge);
356:           if (PCTFS_my_id>dest)
357:           {MPI_Send(vals+segs[edge],stage_n,MPIU_SCALAR,dest,type, MPI_COMM_WORLD);}
358:           else
359:             {
360:               type =  type - PCTFS_my_id + dest;
361:               MPI_Recv(work,stage_n,MPIU_SCALAR,MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
362:               PCTFS_rvec_add(vals+segs[edge], work, stage_n);
363:             }
364:         }
365:       mask <<= 1;
366:     }
367:   mask>>=1;
368:   for (edge=0; edge<level; edge++)
369:     {
370:       stage_n = (segs[level] - segs[level-1-edge]);
371:       if (stage_n && !(PCTFS_my_id & mask))
372:         {
373:           dest = edge_node[level-edge-1];
374:           type = MSGTAG6 + PCTFS_my_id + (PCTFS_num_nodes*edge);
375:           if (PCTFS_my_id<dest)
376:             {MPI_Send(vals+segs[level-1-edge],stage_n,MPIU_SCALAR,dest,type,MPI_COMM_WORLD);}
377:           else
378:             {
379:               type =  type - PCTFS_my_id + dest;
380:               MPI_Recv(vals+segs[level-1-edge],stage_n,MPIU_SCALAR, MPI_ANY_SOURCE,type,MPI_COMM_WORLD,&status);
381:             }
382:         }
383:       mask >>= 1;
384:     }
385:   return(0);
386: }

388: /***********************************comm.c*************************************/
389: PetscErrorCode PCTFS_giop_hc(PetscInt *vals, PetscInt *work, PetscInt n, PetscInt *oprs, PetscInt dim)
390: {
391:   PetscInt            mask, edge;
392:   PetscInt            type, dest;
393:   vfp            fp;
394:   MPI_Status     status;

398:   /* ok ... should have some data, work, and operator(s) */
399:   if (!vals||!work||!oprs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: vals=%D, work=%D, oprs=%D",vals,work,oprs);

401:   /* non-uniform should have at least two entries */
402:   if ((oprs[0] == NON_UNIFORM)&&(n<2)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: non_uniform and n=0,1?");

404:   /* check to make sure comm package has been initialized */
405:   if (!p_init)
406:     {PCTFS_comm_init();}

408:   /* if there's nothing to do return */
409:   if ((PCTFS_num_nodes<2)||(!n)||(dim<=0))
410:     {  return(0);}

412:   /* the error msg says it all!!! */
413:   if (modfl_num_nodes) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: PCTFS_num_nodes not a power of 2!?!");

415:   /* a negative number of items to send ==> fatal */
416:   if (n<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PCTFS_giop_hc() :: n=%D<0?",n);

418:   /* can't do more dimensions then exist */
419:   dim = PetscMin(dim,PCTFS_i_log2_num_nodes);

421:   /* advance to list of n operations for custom */
422:   if ((type=oprs[0])==NON_UNIFORM)
423:     {oprs++;}

425:   if (!(fp = (vfp) PCTFS_ivec_fct_addr(type))){
426:     PetscInfo(0,"PCTFS_giop_hc() :: hope you passed in a rbfp!\n");
427:     fp = (vfp) oprs;
428:   }

430:   for (mask=1,edge=0; edge<dim; edge++,mask<<=1)
431:     {
432:       dest = PCTFS_my_id^mask;
433:       if (PCTFS_my_id > dest)
434:         {MPI_Send(vals,n,MPIU_INT,dest,MSGTAG2+PCTFS_my_id,MPI_COMM_WORLD);}
435:       else
436:         {
437:           MPI_Recv(work,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG2+dest,MPI_COMM_WORLD, &status);
438:           (*fp)(vals, work, n, oprs);
439:         }
440:     }

442:   if (edge==dim)
443:     {mask>>=1;}
444:   else
445:     {while (++edge<dim) {mask<<=1;}}

447:   for (edge=0; edge<dim; edge++,mask>>=1)
448:     {
449:       if (PCTFS_my_id%mask)
450:         {continue;}
451: 
452:       dest = PCTFS_my_id^mask;
453:       if (PCTFS_my_id < dest)
454:         {MPI_Send(vals,n,MPIU_INT,dest,MSGTAG4+PCTFS_my_id,MPI_COMM_WORLD);}
455:       else
456:         {
457:           MPI_Recv(vals,n,MPIU_INT,MPI_ANY_SOURCE,MSGTAG4+dest,MPI_COMM_WORLD,&status);
458:         }
459:     }
460:   return(0);
461: }