Actual source code: xyt.c

  1: #define PETSCKSP_DLL

  3: /*************************************xyt.c************************************
  4: Module Name: xyt
  5: Module Info:

  7: author:  Henry M. Tufo III
  8: e-mail:  hmt@asci.uchicago.edu
  9: contact:
 10: +--------------------------------+--------------------------------+
 11: |MCS Division - Building 221     |Department of Computer Science  |
 12: |Argonne National Laboratory     |Ryerson 152                     |
 13: |9700 S. Cass Avenue             |The University of Chicago       |
 14: |Argonne, IL  60439              |Chicago, IL  60637              |
 15: |(630) 252-5354/5986 ph/fx       |(773) 702-6019/8487 ph/fx       |
 16: +--------------------------------+--------------------------------+

 18: Last Modification: 3.20.01
 19: **************************************xyt.c***********************************/
 20:  #include ../src/ksp/pc/impls/tfs/tfs.h

 22: #define LEFT  -1
 23: #define RIGHT  1
 24: #define BOTH   0

 26: typedef struct xyt_solver_info {
 27:   PetscInt n, m, n_global, m_global;
 28:   PetscInt nnz, max_nnz, msg_buf_sz;
 29:   PetscInt *nsep, *lnsep, *fo, nfo, *stages;
 30:   PetscInt *xcol_sz, *xcol_indices;
 31:   PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
 32:   PetscInt *ycol_sz, *ycol_indices;
 33:   PetscScalar **ycol_vals, *y;
 34:   PetscInt nsolves;
 35:   PetscScalar tot_solve_time;
 36: } xyt_info;

 38: 
 39: typedef struct matvec_info {
 40:   PetscInt n, m, n_global, m_global;
 41:   PetscInt *local2global;
 42:   gs_ADT gs_handle;
 43:   PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
 44:   void *grid_data;
 45: } mv_info;

 47: struct xyt_CDT{
 48:   PetscInt id;
 49:   PetscInt ns;
 50:   PetscInt level;
 51:   xyt_info *info;
 52:   mv_info  *mvi;
 53: };

 55: static PetscInt n_xyt=0;
 56: static PetscInt n_xyt_handles=0;

 58: /* prototypes */
 59: static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
 60: static PetscErrorCode check_handle(xyt_ADT xyt_handle);
 61: static PetscErrorCode det_separators(xyt_ADT xyt_handle);
 62: static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
 63: static PetscInt xyt_generate(xyt_ADT xyt_handle);
 64: static PetscInt do_xyt_factor(xyt_ADT xyt_handle);
 65: static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data);

 67: /**************************************xyt.c***********************************/
 68: xyt_ADT XYT_new(void)
 69: {
 70:   xyt_ADT xyt_handle;

 72:   /* rolling count on n_xyt ... pot. problem here */
 73:   n_xyt_handles++;
 74:   xyt_handle       = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
 75:   xyt_handle->id   = ++n_xyt;
 76:   xyt_handle->info = NULL;
 77:   xyt_handle->mvi  = NULL;

 79:   return(xyt_handle);
 80: }

 82: /**************************************xyt.c***********************************/
 83: PetscInt XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt  handle */
 84:            PetscInt *local2global,  /* global column mapping       */
 85:            PetscInt n,              /* local num rows              */
 86:            PetscInt m,              /* local num cols              */
 87:            void *matvec,       /* b_loc=A_local.x_loc         */
 88:            void *grid_data     /* grid data for matvec        */
 89:            )
 90: {

 92:   comm_init();
 93:   check_handle(xyt_handle);

 95:   /* only 2^k for now and all nodes participating */
 96:   if ((1<<(xyt_handle->level=i_log2_num_nodes))!=num_nodes)
 97:     {SETERRQ2(PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D\n",1<<i_log2_num_nodes,num_nodes);}

 99:   /* space for X info */
100:   xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));

102:   /* set up matvec handles */
103:   xyt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);

105:   /* matrix is assumed to be of full rank */
106:   /* LATER we can reset to indicate rank def. */
107:   xyt_handle->ns=0;

109:   /* determine separators and generate firing order - NB xyt info set here */
110:   det_separators(xyt_handle);

112:   return(do_xyt_factor(xyt_handle));
113: }

115: /**************************************xyt.c***********************************/
116: PetscInt XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
117: {
118:   comm_init();
119:   check_handle(xyt_handle);

121:   /* need to copy b into x? */
122:   if (b)
123:     {rvec_copy(x,b,xyt_handle->mvi->n);}
124:   do_xyt_solve(xyt_handle,x);

126:   return(0);
127: }

129: /**************************************xyt.c***********************************/
130: PetscInt XYT_free(xyt_ADT xyt_handle)
131: {
132:   comm_init();
133:   check_handle(xyt_handle);
134:   n_xyt_handles--;

136:   free(xyt_handle->info->nsep);
137:   free(xyt_handle->info->lnsep);
138:   free(xyt_handle->info->fo);
139:   free(xyt_handle->info->stages);
140:   free(xyt_handle->info->solve_uu);
141:   free(xyt_handle->info->solve_w);
142:   free(xyt_handle->info->x);
143:   free(xyt_handle->info->xcol_vals);
144:   free(xyt_handle->info->xcol_sz);
145:   free(xyt_handle->info->xcol_indices);
146:   free(xyt_handle->info->y);
147:   free(xyt_handle->info->ycol_vals);
148:   free(xyt_handle->info->ycol_sz);
149:   free(xyt_handle->info->ycol_indices);
150:   free(xyt_handle->info);
151:   free(xyt_handle->mvi->local2global);
152:   gs_free(xyt_handle->mvi->gs_handle);
153:   free(xyt_handle->mvi);
154:   free(xyt_handle);

156: 
157:   /* if the check fails we nuke */
158:   /* if NULL pointer passed to free we nuke */
159:   /* if the calls to free fail that's not my problem */
160:   return(0);
161: }

163: /**************************************xyt.c***********************************/
164: PetscInt XYT_stats(xyt_ADT xyt_handle)
165: {
166:   PetscInt  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
167:   PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
168:   PetscInt   vals[9],  work[9];
169:   PetscScalar fvals[3], fwork[3];

172:   comm_init();
173:   check_handle(xyt_handle);

175:   /* if factorization not done there are no stats */
176:   if (!xyt_handle->info||!xyt_handle->mvi)
177:     {
178:       if (!my_id)
179:         {PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n");}
180:       return 1;
181:     }

183:   vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
184:   vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
185:   vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
186:   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);

188:   fvals[0]=fvals[1]=fvals[2]
189:     =xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
190:   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);

192:   if (!my_id)
193:     {
194:       PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_nnz=%D\n",my_id,vals[0]);
195:       PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_nnz=%D\n",my_id,vals[1]);
196:       PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
197:       PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_nnz=%D\n",my_id,vals[2]);
198:       PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
199:       PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
200:       PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_n  =%D\n",my_id,vals[3]);
201:       PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_n  =%D\n",my_id,vals[4]);
202:       PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
203:       PetscPrintf(PETSC_COMM_WORLD,"%D :: tot   xyt_n  =%D\n",my_id,vals[5]);
204:       PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_buf=%D\n",my_id,vals[6]);
205:       PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_buf=%D\n",my_id,vals[7]);
206:       PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
207:       PetscPrintf(PETSC_COMM_WORLD,"%D :: min   xyt_slv=%g\n",my_id,fvals[0]);
208:       PetscPrintf(PETSC_COMM_WORLD,"%D :: max   xyt_slv=%g\n",my_id,fvals[1]);
209:       PetscPrintf(PETSC_COMM_WORLD,"%D :: avg   xyt_slv=%g\n",my_id,fvals[2]/num_nodes);
210:     }

212:   return(0);
213: }


216: /*************************************xyt.c************************************

218: Description: get A_local, local portion of global coarse matrix which 
219: is a row dist. nxm matrix w/ n<m.
220:    o my_ml holds address of ML struct associated w/A_local and coarse grid
221:    o local2global holds global number of column i (i=0,...,m-1)
222:    o local2global holds global number of row    i (i=0,...,n-1)
223:    o mylocmatvec performs A_local . vec_local (note that gs is performed using 
224:    gs_init/gop).

226: mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
227: mylocmatvec (void :: void *data, double *in, double *out)
228: **************************************xyt.c***********************************/
229: static PetscInt do_xyt_factor(xyt_ADT xyt_handle)
230: {
231:   return xyt_generate(xyt_handle);
232: }

234: /**************************************xyt.c***********************************/
235: static PetscInt xyt_generate(xyt_ADT xyt_handle)
236: {
237:   PetscInt i,j,k,idx;
238:   PetscInt dim, col;
239:   PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
240:   PetscInt *segs;
241:   PetscInt op[] = {GL_ADD,0};
242:   PetscInt off, len;
243:   PetscScalar *x_ptr, *y_ptr;
244:   PetscInt *iptr, flag;
245:   PetscInt start=0, end, work;
246:   PetscInt op2[] = {GL_MIN,0};
247:   gs_ADT gs_handle;
248:   PetscInt *nsep, *lnsep, *fo;
249:   PetscInt a_n=xyt_handle->mvi->n;
250:   PetscInt a_m=xyt_handle->mvi->m;
251:   PetscInt *a_local2global=xyt_handle->mvi->local2global;
252:   PetscInt level;
253:   PetscInt n, m;
254:   PetscInt *xcol_sz, *xcol_indices, *stages;
255:   PetscScalar **xcol_vals, *x;
256:   PetscInt *ycol_sz, *ycol_indices;
257:   PetscScalar **ycol_vals, *y;
258:   PetscInt n_global;
259:   PetscInt xt_nnz=0, xt_max_nnz=0;
260:   PetscInt yt_nnz=0, yt_max_nnz=0;
261:   PetscInt xt_zero_nnz  =0;
262:   PetscInt xt_zero_nnz_0=0;
263:   PetscInt yt_zero_nnz  =0;
264:   PetscInt yt_zero_nnz_0=0;
265:   PetscBLASInt i1 = 1,dlen;
266:   PetscScalar dm1 = -1.0;

269:   n=xyt_handle->mvi->n;
270:   nsep=xyt_handle->info->nsep;
271:   lnsep=xyt_handle->info->lnsep;
272:   fo=xyt_handle->info->fo;
273:   end=lnsep[0];
274:   level=xyt_handle->level;
275:   gs_handle=xyt_handle->mvi->gs_handle;

277:   /* is there a null space? */
278:   /* LATER add in ability to detect null space by checking alpha */
279:   for (i=0, j=0; i<=level; i++)
280:     {j+=nsep[i];}

282:   m = j-xyt_handle->ns;
283:   if (m!=j)
284:     {PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);}

286:   PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);

288:   /* get and initialize storage for x local         */
289:   /* note that x local is nxm and stored by columns */
290:   xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
291:   xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
292:   xcol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
293:   for (i=j=0; i<m; i++, j+=2)
294:     {
295:       xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
296:       xcol_vals[i] = NULL;
297:     }
298:   xcol_indices[j]=-1;

300:   /* get and initialize storage for y local         */
301:   /* note that y local is nxm and stored by columns */
302:   ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
303:   ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
304:   ycol_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *));
305:   for (i=j=0; i<m; i++, j+=2)
306:     {
307:       ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
308:       ycol_vals[i] = NULL;
309:     }
310:   ycol_indices[j]=-1;

312:   /* size of separators for each sub-hc working from bottom of tree to top */
313:   /* this looks like nsep[]=segments */
314:   stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
315:   segs   = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
316:   ivec_zero(stages,level+1);
317:   ivec_copy(segs,nsep,level+1);
318:   for (i=0; i<level; i++)
319:     {segs[i+1] += segs[i];}
320:   stages[0] = segs[0];

322:   /* temporary vectors  */
323:   u  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
324:   z  = (PetscScalar *) malloc(n*sizeof(PetscScalar));
325:   v  = (PetscScalar *) malloc(a_m*sizeof(PetscScalar));
326:   uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
327:   w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));

329:   /* extra nnz due to replication of vertices across separators */
330:   for (i=1, j=0; i<=level; i++)
331:     {j+=nsep[i];}

333:   /* storage for sparse x values */
334:   n_global = xyt_handle->info->n_global;
335:   xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
336:   x = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
337:   y = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));

339:   /* LATER - can embed next sep to fire in gs */
340:   /* time to make the donuts - generate X factor */
341:   for (dim=i=j=0;i<m;i++)
342:     {
343:       /* time to move to the next level? */
344:       while (i==segs[dim]){
345:         if (dim==level) SETERRQ(PETSC_ERR_PLIB,"dim about to exceed level\n");
346:         stages[dim++]=i;
347:         end+=lnsep[dim];
348:       }
349:       stages[dim]=i;

351:       /* which column are we firing? */
352:       /* i.e. set v_l */
353:       /* use new seps and do global min across hc to determine which one to fire */
354:       (start<end) ? (col=fo[start]) : (col=INT_MAX);
355:       giop_hc(&col,&work,1,op2,dim);

357:       /* shouldn't need this */
358:       if (col==INT_MAX)
359:         {
360:           PetscInfo(0,"hey ... col==INT_MAX??\n");
361:           continue;
362:         }

364:       /* do I own it? I should */
365:       rvec_zero(v ,a_m);
366:       if (col==fo[start])
367:         {
368:           start++;
369:           idx=ivec_linear_search(col, a_local2global, a_n);
370:           if (idx!=-1)
371:             {v[idx] = 1.0; j++;}
372:           else
373:             {SETERRQ(PETSC_ERR_PLIB,"NOT FOUND!\n");}
374:         }
375:       else
376:         {
377:           idx=ivec_linear_search(col, a_local2global, a_m);
378:           if (idx!=-1)
379:             {v[idx] = 1.0;}
380:         }

382:       /* perform u = A.v_l */
383:       rvec_zero(u,n);
384:       do_matvec(xyt_handle->mvi,v,u);

386:       /* uu =  X^T.u_l (local portion) */
387:       /* technically only need to zero out first i entries */
388:       /* later turn this into an XYT_solve call ? */
389:       rvec_zero(uu,m);
390:       y_ptr=y;
391:       iptr = ycol_indices;
392:       for (k=0; k<i; k++)
393:         {
394:           off = *iptr++;
395:           len = *iptr++;
396:           dlen = PetscBLASIntCast(len);
397:           uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1);
398:           y_ptr+=len;
399:         }

401:       /* uu = X^T.u_l (comm portion) */
402:       ssgl_radd  (uu, w, dim, stages);

404:       /* z = X.uu */
405:       rvec_zero(z,n);
406:       x_ptr=x;
407:       iptr = xcol_indices;
408:       for (k=0; k<i; k++)
409:         {
410:           off = *iptr++;
411:           len = *iptr++;
412:           dlen = PetscBLASIntCast(len);
413:           BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1);
414:           x_ptr+=len;
415:         }

417:       /* compute v_l = v_l - z */
418:       rvec_zero(v+a_n,a_m-a_n);
419:       dlen = PetscBLASIntCast(n);
420:       BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1);

422:       /* compute u_l = A.v_l */
423:       if (a_n!=a_m)
424:         {gs_gop_hc(gs_handle,v,"+\0",dim);}
425:       rvec_zero(u,n);
426:      do_matvec(xyt_handle->mvi,v,u);

428:       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
429:      dlen = PetscBLASIntCast(n);
430:       alpha = BLASdot_(&dlen,u,&i1,u,&i1);
431:       /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
432:       grop_hc(&alpha, &alpha_w, 1, op, dim);

434:       alpha = (PetscScalar) sqrt((double)alpha);

436:       /* check for small alpha                             */
437:       /* LATER use this to detect and determine null space */
438:       if (fabs(alpha)<1.0e-14)
439:         {SETERRQ1(PETSC_ERR_PLIB,"bad alpha! %g\n",alpha);}

441:       /* compute v_l = v_l/sqrt(alpha) */
442:       rvec_scale(v,1.0/alpha,n);
443:       rvec_scale(u,1.0/alpha,n);

445:       /* add newly generated column, v_l, to X */
446:       flag = 1;
447:       off=len=0;
448:       for (k=0; k<n; k++)
449:         {
450:           if (v[k]!=0.0)
451:             {
452:               len=k;
453:               if (flag)
454:                 {off=k; flag=0;}
455:             }
456:         }

458:       len -= (off-1);

460:       if (len>0)
461:         {
462:           if ((xt_nnz+len)>xt_max_nnz)
463:             {
464:               PetscInfo(0,"increasing space for X by 2x!\n");
465:               xt_max_nnz *= 2;
466:               x_ptr = (PetscScalar *) malloc(xt_max_nnz*sizeof(PetscScalar));
467:               rvec_copy(x_ptr,x,xt_nnz);
468:               free(x);
469:               x = x_ptr;
470:               x_ptr+=xt_nnz;
471:             }
472:           xt_nnz += len;
473:           rvec_copy(x_ptr,v+off,len);

475:           /* keep track of number of zeros */
476:           if (dim)
477:             {
478:               for (k=0; k<len; k++)
479:                 {
480:                   if (x_ptr[k]==0.0)
481:                     {xt_zero_nnz++;}
482:                 }
483:             }
484:           else
485:             {
486:               for (k=0; k<len; k++)
487:                 {
488:                   if (x_ptr[k]==0.0)
489:                     {xt_zero_nnz_0++;}
490:                 }
491:             }
492:           xcol_indices[2*i] = off;
493:           xcol_sz[i] = xcol_indices[2*i+1] = len;
494:           xcol_vals[i] = x_ptr;
495:         }
496:       else
497:         {
498:           xcol_indices[2*i] = 0;
499:           xcol_sz[i] = xcol_indices[2*i+1] = 0;
500:           xcol_vals[i] = x_ptr;
501:         }


504:       /* add newly generated column, u_l, to Y */
505:       flag = 1;
506:       off=len=0;
507:       for (k=0; k<n; k++)
508:         {
509:           if (u[k]!=0.0)
510:             {
511:               len=k;
512:               if (flag)
513:                 {off=k; flag=0;}
514:             }
515:         }

517:       len -= (off-1);

519:       if (len>0)
520:         {
521:           if ((yt_nnz+len)>yt_max_nnz)
522:             {
523:               PetscInfo(0,"increasing space for Y by 2x!\n");
524:               yt_max_nnz *= 2;
525:               y_ptr = (PetscScalar *) malloc(yt_max_nnz*sizeof(PetscScalar));
526:               rvec_copy(y_ptr,y,yt_nnz);
527:               free(y);
528:               y = y_ptr;
529:               y_ptr+=yt_nnz;
530:             }
531:           yt_nnz += len;
532:           rvec_copy(y_ptr,u+off,len);

534:           /* keep track of number of zeros */
535:           if (dim)
536:             {
537:               for (k=0; k<len; k++)
538:                 {
539:                   if (y_ptr[k]==0.0)
540:                     {yt_zero_nnz++;}
541:                 }
542:             }
543:           else
544:             {
545:               for (k=0; k<len; k++)
546:                 {
547:                   if (y_ptr[k]==0.0)
548:                     {yt_zero_nnz_0++;}
549:                 }
550:             }
551:           ycol_indices[2*i] = off;
552:           ycol_sz[i] = ycol_indices[2*i+1] = len;
553:           ycol_vals[i] = y_ptr;
554:         }
555:       else
556:         {
557:           ycol_indices[2*i] = 0;
558:           ycol_sz[i] = ycol_indices[2*i+1] = 0;
559:           ycol_vals[i] = y_ptr;
560:         }
561:     }

563:   /* close off stages for execution phase */
564:   while (dim!=level)
565:     {
566:       stages[dim++]=i;
567:       PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);
568:     }
569:   stages[dim]=i;

571:   xyt_handle->info->n=xyt_handle->mvi->n;
572:   xyt_handle->info->m=m;
573:   xyt_handle->info->nnz=xt_nnz + yt_nnz;
574:   xyt_handle->info->max_nnz=xt_max_nnz + yt_max_nnz;
575:   xyt_handle->info->msg_buf_sz=stages[level]-stages[0];
576:   xyt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar));
577:   xyt_handle->info->solve_w  = (PetscScalar *) malloc(m*sizeof(PetscScalar));
578:   xyt_handle->info->x=x;
579:   xyt_handle->info->xcol_vals=xcol_vals;
580:   xyt_handle->info->xcol_sz=xcol_sz;
581:   xyt_handle->info->xcol_indices=xcol_indices;
582:   xyt_handle->info->stages=stages;
583:   xyt_handle->info->y=y;
584:   xyt_handle->info->ycol_vals=ycol_vals;
585:   xyt_handle->info->ycol_sz=ycol_sz;
586:   xyt_handle->info->ycol_indices=ycol_indices;

588:   free(segs);
589:   free(u);
590:   free(v);
591:   free(uu);
592:   free(z);
593:   free(w);

595:   return(0);
596: }

598: /**************************************xyt.c***********************************/
599: static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle,  PetscScalar *uc)
600: {
601:   PetscInt off, len, *iptr;
602:   PetscInt level       =xyt_handle->level;
603:   PetscInt n           =xyt_handle->info->n;
604:   PetscInt m           =xyt_handle->info->m;
605:   PetscInt *stages     =xyt_handle->info->stages;
606:   PetscInt *xcol_indices=xyt_handle->info->xcol_indices;
607:   PetscInt *ycol_indices=xyt_handle->info->ycol_indices;
608:    PetscScalar *x_ptr, *y_ptr, *uu_ptr;
609:   PetscScalar *solve_uu=xyt_handle->info->solve_uu;
610:   PetscScalar *solve_w =xyt_handle->info->solve_w;
611:   PetscScalar *x       =xyt_handle->info->x;
612:   PetscScalar *y       =xyt_handle->info->y;
613:   PetscBLASInt i1 = 1,dlen;

616:   uu_ptr=solve_uu;
617:   rvec_zero(uu_ptr,m);

619:   /* x  = X.Y^T.b */
620:   /* uu = Y^T.b */
621:   for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len)
622:     {
623:       off=*iptr++; len=*iptr++;
624:       dlen = PetscBLASIntCast(len);
625:       *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1);
626:     }

628:   /* comunication of beta */
629:   uu_ptr=solve_uu;
630:   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}

632:   rvec_zero(uc,n);

634:   /* x = X.uu */
635:   for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len)
636:     {
637:       off=*iptr++; len=*iptr++;
638:       dlen = PetscBLASIntCast(len);
639:       BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1);
640:     }
641:   return(0);
642: }

644: /**************************************xyt.c***********************************/
645: static PetscErrorCode check_handle(xyt_ADT xyt_handle)
646: {
647:   PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};

650:   if (xyt_handle==NULL)
651:     {SETERRQ1(PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle);}

653:   vals[0]=vals[1]=xyt_handle->id;
654:   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
655:   if ((vals[0]!=vals[1])||(xyt_handle->id<=0))
656:     {SETERRQ3(PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D\n", vals[0],vals[1], xyt_handle->id);}
657:   return(0);
658: }

660: /**************************************xyt.c***********************************/
661: static PetscErrorCode det_separators(xyt_ADT xyt_handle)
662: {
663:   PetscInt i, ct, id;
664:   PetscInt mask, edge, *iptr;
665:   PetscInt *dir, *used;
666:   PetscInt sum[4], w[4];
667:   PetscScalar rsum[4], rw[4];
668:   PetscInt op[] = {GL_ADD,0};
669:   PetscScalar *lhs, *rhs;
670:   PetscInt *nsep, *lnsep, *fo, nfo=0;
671:   gs_ADT gs_handle=xyt_handle->mvi->gs_handle;
672:   PetscInt *local2global=xyt_handle->mvi->local2global;
673:   PetscInt  n=xyt_handle->mvi->n;
674:   PetscInt  m=xyt_handle->mvi->m;
675:   PetscInt level=xyt_handle->level;
676:   PetscInt shared=FALSE;

680:   dir  = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
681:   nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
682:   lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
683:   fo   = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
684:   used = (PetscInt*)malloc(sizeof(PetscInt)*n);

686:   ivec_zero(dir  ,level+1);
687:   ivec_zero(nsep ,level+1);
688:   ivec_zero(lnsep,level+1);
689:   ivec_set (fo   ,-1,n+1);
690:   ivec_zero(used,n);

692:   lhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
693:   rhs  = (PetscScalar*)malloc(sizeof(PetscScalar)*m);

695:   /* determine the # of unique dof */
696:   rvec_zero(lhs,m);
697:   rvec_set(lhs,1.0,n);
698:   gs_gop_hc(gs_handle,lhs,"+\0",level);
699:   PetscInfo(0,"done first gs_gop_hc\n");
700:   rvec_zero(rsum,2);
701:   for (ct=i=0;i<n;i++)
702:     {
703:       if (lhs[i]!=0.0)
704:         {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}

706:       if (lhs[i]!=1.0)
707:         {
708:           shared=TRUE;
709:         }
710:     }

712:   grop_hc(rsum,rw,2,op,level);
713:   rsum[0]+=0.1;
714:   rsum[1]+=0.1;

716:   xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
717:   xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];

719:   /* determine separator sets top down */
720:   if (shared)
721:     {
722:       /* solution is to do as in the symmetric shared case but then */
723:       /* pick the sub-hc with the most free dofs and do a mat-vec   */
724:       /* and pick up the responses on the other sub-hc from the     */
725:       /* initial separator set obtained from the symm. shared case  */
726:       SETERRQ(PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n");
727:       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
728:         {
729:           /* set rsh of hc, fire, and collect lhs responses */
730:           (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
731:           gs_gop_hc(gs_handle,lhs,"+\0",edge);
732: 
733:           /* set lsh of hc, fire, and collect rhs responses */
734:           (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
735:           gs_gop_hc(gs_handle,rhs,"+\0",edge);
736: 
737:           for (i=0;i<n;i++)
738:             {
739:               if (id< mask)
740:                 {
741:                   if (lhs[i]!=0.0)
742:                     {lhs[i]=1.0;}
743:                 }
744:               if (id>=mask)
745:                 {
746:                   if (rhs[i]!=0.0)
747:                     {rhs[i]=1.0;}
748:                 }
749:             }

751:           if (id< mask)
752:             {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
753:           else
754:             {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}

756:           /* count number of dofs I own that have signal and not in sep set */
757:           rvec_zero(rsum,4);
758:           for (ivec_zero(sum,4),ct=i=0;i<n;i++)
759:             {
760:               if (!used[i])
761:                 {
762:                   /* number of unmarked dofs on node */
763:                   ct++;
764:                   /* number of dofs to be marked on lhs hc */
765:                   if (id< mask)
766:                     {
767:                       if (lhs[i]!=0.0)
768:                         {sum[0]++; rsum[0]+=1.0/lhs[i];}
769:                     }
770:                   /* number of dofs to be marked on rhs hc */
771:                   if (id>=mask)
772:                     {
773:                       if (rhs[i]!=0.0)
774:                         {sum[1]++; rsum[1]+=1.0/rhs[i];}
775:                     }
776:                 }
777:             }

779:           /* go for load balance - choose half with most unmarked dofs, bias LHS */
780:           (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
781:           (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
782:           giop_hc(sum,w,4,op,edge);
783:           grop_hc(rsum,rw,4,op,edge);
784:           rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;

786:           if (id<mask)
787:             {
788:               /* mark dofs I own that have signal and not in sep set */
789:               for (ct=i=0;i<n;i++)
790:                 {
791:                   if ((!used[i])&&(lhs[i]!=0.0))
792:                     {
793:                       ct++; nfo++;

795:                       if (nfo>n)
796:                         {SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}

798:                       *--iptr = local2global[i];
799:                       used[i]=edge;
800:                     }
801:                 }
802:               if (ct>1) {ivec_sort(iptr,ct);}

804:               lnsep[edge]=ct;
805:               nsep[edge]=(PetscInt) rsum[0];
806:               dir [edge]=LEFT;
807:             }

809:           if (id>=mask)
810:             {
811:               /* mark dofs I own that have signal and not in sep set */
812:               for (ct=i=0;i<n;i++)
813:                 {
814:                   if ((!used[i])&&(rhs[i]!=0.0))
815:                     {
816:                       ct++; nfo++;

818:                       if (nfo>n)
819:                         {SETERRQ(PETSC_ERR_PLIB,"nfo about to exceed n\n");}

821:                       *--iptr = local2global[i];
822:                       used[i]=edge;
823:                     }
824:                 }
825:               if (ct>1) {ivec_sort(iptr,ct);}

827:               lnsep[edge]=ct;
828:               nsep[edge]= (PetscInt) rsum[1];
829:               dir [edge]=RIGHT;
830:             }

832:           /* LATER or we can recur on these to order seps at this level */
833:           /* do we need full set of separators for this?                */

835:           /* fold rhs hc into lower */
836:           if (id>=mask)
837:             {id-=mask;}
838:         }
839:     }
840:   else
841:     {
842:       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
843:         {
844:           /* set rsh of hc, fire, and collect lhs responses */
845:           (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
846:           gs_gop_hc(gs_handle,lhs,"+\0",edge);

848:           /* set lsh of hc, fire, and collect rhs responses */
849:           (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
850:           gs_gop_hc(gs_handle,rhs,"+\0",edge);

852:           /* count number of dofs I own that have signal and not in sep set */
853:           for (ivec_zero(sum,4),ct=i=0;i<n;i++)
854:             {
855:               if (!used[i])
856:                 {
857:                   /* number of unmarked dofs on node */
858:                   ct++;
859:                   /* number of dofs to be marked on lhs hc */
860:                   if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
861:                   /* number of dofs to be marked on rhs hc */
862:                   if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
863:                 }
864:             }

866:           /* for the non-symmetric case we need separators of width 2 */
867:           /* so take both sides */
868:           (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
869:           giop_hc(sum,w,4,op,edge);

871:           ct=0;
872:           if (id<mask)
873:             {
874:               /* mark dofs I own that have signal and not in sep set */
875:               for (i=0;i<n;i++)
876:                 {
877:                   if ((!used[i])&&(lhs[i]!=0.0))
878:                     {
879:                       ct++; nfo++;
880:                       *--iptr = local2global[i];
881:                       used[i]=edge;
882:                     }
883:                 }
884:               /* LSH hc summation of ct should be sum[0] */
885:             }
886:           else
887:             {
888:               /* mark dofs I own that have signal and not in sep set */
889:               for (i=0;i<n;i++)
890:                 {
891:                   if ((!used[i])&&(rhs[i]!=0.0))
892:                     {
893:                       ct++; nfo++;
894:                       *--iptr = local2global[i];
895:                       used[i]=edge;
896:                     }
897:                 }
898:               /* RSH hc summation of ct should be sum[1] */
899:             }

901:           if (ct>1) {ivec_sort(iptr,ct);}
902:           lnsep[edge]=ct;
903:           nsep[edge]=sum[0]+sum[1];
904:           dir [edge]=BOTH;

906:           /* LATER or we can recur on these to order seps at this level */
907:           /* do we need full set of separators for this?                */

909:           /* fold rhs hc into lower */
910:           if (id>=mask)
911:             {id-=mask;}
912:         }
913:     }

915:   /* level 0 is on processor case - so mark the remainder */
916:   for (ct=i=0;i<n;i++)
917:     {
918:       if (!used[i])
919:         {
920:           ct++; nfo++;
921:           *--iptr = local2global[i];
922:           used[i]=edge;
923:         }
924:     }
925:   if (ct>1) {ivec_sort(iptr,ct);}
926:   lnsep[edge]=ct;
927:   nsep [edge]=ct;
928:   dir  [edge]=BOTH;

930:   xyt_handle->info->nsep=nsep;
931:   xyt_handle->info->lnsep=lnsep;
932:   xyt_handle->info->fo=fo;
933:   xyt_handle->info->nfo=nfo;

935:   free(dir);
936:   free(lhs);
937:   free(rhs);
938:   free(used);
939:   return(0);
940: }

942: /**************************************xyt.c***********************************/
943: static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data)
944: {
945:   mv_info *mvi;


948:   mvi = (mv_info*)malloc(sizeof(mv_info));
949:   mvi->n=n;
950:   mvi->m=m;
951:   mvi->n_global=-1;
952:   mvi->m_global=-1;
953:   mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt));
954:   ivec_copy(mvi->local2global,local2global,m);
955:   mvi->local2global[m] = INT_MAX;
956:   mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec;
957:   mvi->grid_data=grid_data;

959:   /* set xyt communication handle to perform restricted matvec */
960:   mvi->gs_handle = gs_init(local2global, m, num_nodes);

962:   return(mvi);
963: }

965: /**************************************xyt.c***********************************/
966: static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
967: {
969:   A->matvec((mv_info*)A->grid_data,v,u);
970:   return(0);
971: }