Actual source code: xxt.c

  1: /*$Id: xxt.c,v 1.5 2001/04/12 15:23:12 bsmith Exp $*/
  2: /*************************************xxt.c************************************
  3: Module Name: xxt
  4: Module Info:

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

 17: Last Modification: 3.20.01
 18: **************************************xxt.c***********************************/


 21: /*************************************xxt.c************************************
 22: NOTES ON USAGE: 

 24: **************************************xxt.c***********************************/
 25: #include <stdio.h>
 26: #include <stdlib.h>
 27: #include <limits.h>
 28: #include <float.h>
 29: #include <math.h>
 30: #if   defined NXSRC
 31: #include <nx.h>
 32: #elif defined MPISRC
 33: #include <mpi.h>
 34: #endif

 36:  #include const.h
 37:  #include types.h
 38:  #include comm.h
 39:  #include error.h
 40:  #include ivec.h
 41: #include "bss_malloc.h"
 42:  #include queue.h
 43:  #include gs.h
 44: #ifdef MLSRC
 45: #include "ml_include.h"
 46: #endif
 47:  #include blas.h
 48:  #include xxt.h

 50: #define LEFT  -1
 51: #define RIGHT  1
 52: #define BOTH   0
 53: #define MAX_FORTRAN_HANDLES  10

 55: typedef struct xxt_solver_info {
 56:   int n, m, n_global, m_global;
 57:   int nnz, max_nnz, msg_buf_sz;
 58:   int *nsep, *lnsep, *fo, nfo, *stages;
 59:   int *col_sz, *col_indices;
 60:   REAL **col_vals, *x, *solve_uu, *solve_w;
 61:   int nsolves;
 62:   REAL tot_solve_time;
 63: } xxt_info;

 65: typedef struct matvec_info {
 66:   int n, m, n_global, m_global;
 67:   int *local2global;
 68:   gs_ADT gs_handle;
 69:   int (*matvec)(struct matvec_info*,REAL*,REAL*);
 70:   void *grid_data;
 71: } mv_info;

 73: struct xxt_CDT{
 74:   int id;
 75:   int ns;
 76:   int level;
 77:   xxt_info *info;
 78:   mv_info  *mvi;
 79: };

 81: static int n_xxt=0;
 82: static int n_xxt_handles=0;

 84: /* prototypes */
 85: static void do_xxt_solve(xxt_ADT xxt_handle, REAL *rhs);
 86: static void check_init(void);
 87: static void check_handle(xxt_ADT xxt_handle);
 88: static void det_separators(xxt_ADT xxt_handle);
 89: static void do_matvec(mv_info *A, REAL *v, REAL *u);
 90: static int xxt_generate(xxt_ADT xxt_handle);
 91: static int do_xxt_factor(xxt_ADT xxt_handle);
 92: static mv_info *set_mvi(int *local2global, int n, int m, void *matvec, void *grid_data);
 93: #ifdef MLSRC
 94: void ML_XXT_solve(xxt_ADT xxt_handle, int lx, double *x, int lb, double *b);
 95: int  ML_XXT_factor(xxt_ADT xxt_handle, int *local2global, int n, int m,
 96:                    void *matvec, void *grid_data, int grid_tag, ML *my_ml);
 97: #endif


100: /*************************************xxt.c************************************
101: Function: XXT_new()

103: Input :
104: Output:
105: Return:
106: Description:
107: **************************************xxt.c***********************************/
108: xxt_ADT 
109: XXT_new(void)
110: {
111:   xxt_ADT xxt_handle;


114: #ifdef DEBUG
115:   error_msg_warning("XXT_new() :: start %d\n",n_xxt_handles);
116: #endif

118:   /* rolling count on n_xxt ... pot. problem here */
119:   n_xxt_handles++;
120:   xxt_handle       = (xxt_ADT)bss_malloc(sizeof(struct xxt_CDT));
121:   xxt_handle->id   = ++n_xxt;
122:   xxt_handle->info = NULL; xxt_handle->mvi  = NULL;

124: #ifdef DEBUG
125:   error_msg_warning("XXT_new() :: end   %d\n",n_xxt_handles);
126: #endif

128:   return(xxt_handle);
129: }


132: /*************************************xxt.c************************************
133: Function: XXT_factor()

135: Input :
136: Output:
137: Return:
138: Description:
139: **************************************xxt.c***********************************/
140: int
141: XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt  handle */
142:            int *local2global,  /* global column mapping       */
143:            int n,              /* local num rows              */
144:            int m,              /* local num cols              */
145:            void *matvec,       /* b_loc=A_local.x_loc         */
146:            void *grid_data     /* grid data for matvec        */
147:            )
148: {
149: #ifdef DEBUG
150:   int flag;


153:   error_msg_warning("XXT_factor() :: start %d\n",n_xxt_handles);
154: #endif

156:   check_init();
157:   check_handle(xxt_handle);

159:   /* only 2^k for now and all nodes participating */
160:   if ((1<<(xxt_handle->level=i_log2_num_nodes))!=num_nodes)
161:     {error_msg_fatal("only 2^k for now and MPI_COMM_WORLD!!! %d != %d\n",1<<i_log2_num_nodes,num_nodes);}

163:   /* space for X info */
164:   xxt_handle->info = (xxt_info*)bss_malloc(sizeof(xxt_info));

166:   /* set up matvec handles */
167:   xxt_handle->mvi  = set_mvi(local2global, n, m, matvec, grid_data);

169:   /* matrix is assumed to be of full rank */
170:   /* LATER we can reset to indicate rank def. */
171:   xxt_handle->ns=0;

173:   /* determine separators and generate firing order - NB xxt info set here */
174:   det_separators(xxt_handle);

176: #ifdef DEBUG
177:   flag = do_xxt_factor(xxt_handle);
178:   error_msg_warning("XXT_factor() :: end   %d (flag=%d)\n",n_xxt_handles,flag);
179:   return(flag);
180: #else
181:   return(do_xxt_factor(xxt_handle));
182: #endif
183: }


186: /*************************************xxt.c************************************
187: Function: XXT_solve

189: Input :
190: Output:
191: Return:
192: Description:
193: **************************************xxt.c***********************************/
194: int
195: XXT_solve(xxt_ADT xxt_handle, double *x, double *b)
196: {
197: #if defined(NXSRC) && defined(TIMING)
198:   double dclock(),    time=0.0;
199: #elif defined(MPISRC) && defined(TIMING)
200:   double MPI_Wtime(), time=0.0;
201: #endif
202: #ifdef INFO
203:   REAL vals[3], work[3];
204:   int op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
205: #endif


208: #ifdef DEBUG
209:   error_msg_warning("XXT_solve() :: start %d\n",n_xxt_handles);
210: #endif

212:   check_init();
213:   check_handle(xxt_handle);

215: #if defined(NXSRC) && defined(TIMING)
216:   time = dclock();
217: #elif defined(MPISRC) && defined(TIMING)
218:   time = MPI_Wtime();
219: #endif

221:   /* need to copy b into x? */
222:   if (b)
223:     {rvec_copy(x,b,xxt_handle->mvi->n);}
224:   do_xxt_solve(xxt_handle,x);

226: #if defined(NXSRC) && defined(TIMING)
227:   time = dclock() - time;
228: #elif defined(MPISRC) && defined(TIMING)
229:   time = MPI_Wtime() - time;
230: #endif

232: #ifdef TIMING
233:   xxt_handle->info->nsolves++;
234:   xxt_handle->info->tot_solve_time+=time;

236: #ifdef INFO
237:   vals[0]=vals[1]=vals[2]= (REAL) time;
238:   grop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
239:   if (!my_id)
240:     {
241:       printf("%d :: min   xxt_slv=%g\n",my_id,vals[0]);
242:       printf("%d :: max   xxt_slv=%g\n",my_id,vals[1]);
243:       printf("%d :: avg   xxt_slv=%g\n",my_id,vals[2]/num_nodes);
244:     }
245: #endif
246: #endif

248: #ifdef DEBUG
249:   error_msg_warning("XXT_solve() :: end   %d\n",n_xxt_handles);
250: #endif

252:   return(0);
253: }


256: /*************************************xxt.c************************************
257: Function: XXT_free()

259: Input :
260: Output:
261: Return:
262: Description:
263: **************************************xxt.c***********************************/
264: int
265: XXT_free(xxt_ADT xxt_handle)
266: {
267: #ifdef DEBUG
268:   error_msg_warning("XXT_free() :: start %d\n",n_xxt_handles);
269: #endif

271:   check_init();
272:   check_handle(xxt_handle);
273:   n_xxt_handles--;

275:   bss_free(xxt_handle->info->nsep);
276:   bss_free(xxt_handle->info->lnsep);
277:   bss_free(xxt_handle->info->fo);
278:   bss_free(xxt_handle->info->stages);
279:   bss_free(xxt_handle->info->solve_uu);
280:   bss_free(xxt_handle->info->solve_w);
281:   bss_free(xxt_handle->info->x);
282:   bss_free(xxt_handle->info->col_vals);
283:   bss_free(xxt_handle->info->col_sz);
284:   bss_free(xxt_handle->info->col_indices);
285:   bss_free(xxt_handle->info);
286:   bss_free(xxt_handle->mvi->local2global);
287:    gs_free(xxt_handle->mvi->gs_handle);
288:   bss_free(xxt_handle->mvi);
289:   bss_free(xxt_handle);

291: 
292: #ifdef DEBUG
293:   error_msg_warning("perm frees = %d\n",perm_frees());
294:   error_msg_warning("perm calls = %d\n",perm_calls());
295:   error_msg_warning("bss frees  = %d\n",bss_frees());
296:   error_msg_warning("bss calls  = %d\n",bss_calls());
297:   error_msg_warning("XXT_free() :: end   %d\n",n_xxt_handles);
298: #endif

300:   /* if the check fails we nuke */
301:   /* if NULL pointer passed to bss_free we nuke */
302:   /* if the calls to free fail that's not my problem */
303:   return(0);
304: }


307: #ifdef MLSRC
308: /*************************************xxt.c************************************
309: Function: ML_XXT_factor()

311: Input :
312: Output:
313: Return:
314: Description:

316: ML requires that the solver call be checked in
317: **************************************xxt.c***********************************/
318: int
319: ML_XXT_factor(xxt_ADT xxt_handle,  /* prev. allocated xxt  handle */
320:                 int *local2global, /* global column mapping       */
321:                 int n,             /* local num rows              */
322:                 int m,             /* local num cols              */
323:                 void *matvec,      /* b_loc=A_local.x_loc         */
324:                 void *grid_data,   /* grid data for matvec        */
325:                 int grid_tag,      /* grid tag for ML_Set_CSolve  */
326:                 ML *my_ml          /* ML handle                   */
327:                 )
328: {
329: #ifdef DEBUG
330:   int flag;
331: #endif


334: #ifdef DEBUG
335:   error_msg_warning("ML_XXT_factor() :: start %d\n",n_xxt_handles);
336: #endif

338:   check_init();
339:   check_handle(xxt_handle);
340:   if (my_ml->comm->ML_mypid!=my_id)
341:     {error_msg_fatal("ML_XXT_factor bad my_id %d\t%d\n",
342:                      my_ml->comm->ML_mypid,my_id);}
343:   if (my_ml->comm->ML_nprocs!=num_nodes)
344:     {error_msg_fatal("ML_XXT_factor bad np %d\t%d\n",
345:                      my_ml->comm->ML_nprocs,num_nodes);}

347:   my_ml->SingleLevel[grid_tag].csolve->func->external = ML_XXT_solve;
348:   my_ml->SingleLevel[grid_tag].csolve->func->ML_id = ML_EXTERNAL;
349:   my_ml->SingleLevel[grid_tag].csolve->data = xxt_handle;

351:   /* done ML specific stuff ... back to reg sched pgm */
352: #ifdef DEBUG
353:   flag = XXT_factor(xxt_handle, local2global, n, m, matvec, grid_data);
354:   error_msg_warning("ML_XXT_factor() :: end   %d (flag=%d)\n",n_xxt_handles,flag);
355:   return(flag);
356: #else
357:   return(XXT_factor(xxt_handle, local2global, n, m, matvec, grid_data));
358: #endif
359: }


362: /*************************************xxt.c************************************
363: Function: ML_XXT_solve

365: Input :
366: Output:
367: Return:
368: Description:
369: **************************************xxt.c***********************************/
370: void 
371: ML_XXT_solve(xxt_ADT xxt_handle, int lx, double *sol, int lb, double *rhs)
372: {
373:   XXT_solve(xxt_handle, sol, rhs);
374: }
375: #endif


378: /*************************************xxt.c************************************
379: Function: 

381: Input : 
382: Output: 
383: Return: 
384: Description:  
385: **************************************xxt.c***********************************/
386: int
387: XXT_stats(xxt_ADT xxt_handle)
388: {
389:   int  op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
390:   int fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
391:   int   vals[9],  work[9];
392:   REAL fvals[3], fwork[3];


395: #ifdef DEBUG
396:   error_msg_warning("xxt_stats() :: begin\n");
397: #endif

399:   check_init();
400:   check_handle(xxt_handle);

402:   /* if factorization not done there are no stats */
403:   if (!xxt_handle->info||!xxt_handle->mvi)
404:     {
405:       if (!my_id)
406:         {printf("XXT_stats() :: no stats available!\n");}
407:       return 1;
408:     }

410:   vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz;
411:   vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n;
412:   vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz;
413:   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);

415:   fvals[0]=fvals[1]=fvals[2]
416:     =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++;
417:   grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);

419:   if (!my_id)
420:     {
421:       printf("%d :: min   xxt_nnz=%d\n",my_id,vals[0]);
422:       printf("%d :: max   xxt_nnz=%d\n",my_id,vals[1]);
423:       printf("%d :: avg   xxt_nnz=%g\n",my_id,1.0*vals[2]/num_nodes);
424:       printf("%d :: tot   xxt_nnz=%d\n",my_id,vals[2]);
425:       printf("%d :: xxt   C(2d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.5)));
426:       printf("%d :: xxt   C(3d)  =%g\n",my_id,vals[2]/(pow(1.0*vals[5],1.6667)));
427:       printf("%d :: min   xxt_n  =%d\n",my_id,vals[3]);
428:       printf("%d :: max   xxt_n  =%d\n",my_id,vals[4]);
429:       printf("%d :: avg   xxt_n  =%g\n",my_id,1.0*vals[5]/num_nodes);
430:       printf("%d :: tot   xxt_n  =%d\n",my_id,vals[5]);
431:       printf("%d :: min   xxt_buf=%d\n",my_id,vals[6]);
432:       printf("%d :: max   xxt_buf=%d\n",my_id,vals[7]);
433:       printf("%d :: avg   xxt_buf=%g\n",my_id,1.0*vals[8]/num_nodes);
434:       printf("%d :: min   xxt_slv=%g\n",my_id,fvals[0]);
435:       printf("%d :: max   xxt_slv=%g\n",my_id,fvals[1]);
436:       printf("%d :: avg   xxt_slv=%g\n",my_id,fvals[2]/num_nodes);
437:     }

439: #ifdef DEBUG
440:   error_msg_warning("xxt_stats() :: end\n");
441: #endif

443:   return(0);
444: }


447: /*************************************xxt.c************************************
448: Function: do_xxt_factor

450: Input : 
451: Output: 
452: Return: 
453: Description: get A_local, local portion of global coarse matrix which 
454: is a row dist. nxm matrix w/ n<m.
455:    o my_ml holds address of ML struct associated w/A_local and coarse grid
456:    o local2global holds global number of column i (i=0,...,m-1)
457:    o local2global holds global number of row    i (i=0,...,n-1)
458:    o mylocmatvec performs A_local . vec_local (note that gs is performed using 
459:    gs_init/gop).

461: mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
462: mylocmatvec (void :: void *data, double *in, double *out)
463: **************************************xxt.c***********************************/
464: static
465: int
466: do_xxt_factor(xxt_ADT xxt_handle)
467: {
468:   int flag;


471: #ifdef DEBUG
472:   error_msg_warning("do_xxt_factor() :: begin\n");
473: #endif

475:   flag=xxt_generate(xxt_handle);

477: #ifdef INFO
478:   XXT_stats(xxt_handle);
479:   bss_stats();
480:   perm_stats();
481: #endif

483: #ifdef DEBUG
484:   error_msg_warning("do_xxt_factor() :: end\n");
485: #endif

487:   return(flag);
488: }


491: /*************************************xxt.c************************************
492: Function: 

494: Input : 
495: Output: 
496: Return: 
497: Description:  
498: **************************************xxt.c***********************************/
499: static
500: int
501: xxt_generate(xxt_ADT xxt_handle)
502: {
503:   int i,j,k,idex;
504:   int dim, col;
505:   REAL *u, *uu, *v, *z, *w, alpha, alpha_w;
506:   int *segs;
507:   int op[] = {GL_ADD,0};
508:   int off, len;
509:   REAL *x_ptr;
510:   int *iptr, flag;
511:   int start=0, end, work;
512:   int op2[] = {GL_MIN,0};
513:   gs_ADT gs_handle;
514:   int *nsep, *lnsep, *fo;
515:   int a_n=xxt_handle->mvi->n;
516:   int a_m=xxt_handle->mvi->m;
517:   int *a_local2global=xxt_handle->mvi->local2global;
518:   int level;
519:   int xxt_nnz=0, xxt_max_nnz=0;
520:   int n, m;
521:   int *col_sz, *col_indices, *stages;
522:   REAL **col_vals, *x;
523:   int n_global;
524:   int xxt_zero_nnz=0;
525:   int xxt_zero_nnz_0=0;


528: #ifdef DEBUG
529:   error_msg_warning("xxt_generate() :: begin\n");
530: #endif

532:   n=xxt_handle->mvi->n;
533:   nsep=xxt_handle->info->nsep;
534:   lnsep=xxt_handle->info->lnsep;
535:   fo=xxt_handle->info->fo;
536:   end=lnsep[0];
537:   level=xxt_handle->level;
538:   gs_handle=xxt_handle->mvi->gs_handle;

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

545:   m = j-xxt_handle->ns;
546:   if (m!=j)
547:     {printf("xxt_generate() :: null space exists %d %d %d\n",m,j,xxt_handle->ns);}

549:   /* get and initialize storage for x local         */
550:   /* note that x local is nxm and stored by columns */
551:   col_sz = (int *) bss_malloc(m*INT_LEN);
552:   col_indices = (int *) bss_malloc((2*m+1)*sizeof(int));
553:   col_vals = (REAL **) bss_malloc(m*sizeof(REAL *));
554:   for (i=j=0; i<m; i++, j+=2)
555:     {
556:       col_indices[j]=col_indices[j+1]=col_sz[i]=-1;
557:       col_vals[i] = NULL;
558:     }
559:   col_indices[j]=-1;

561:   /* size of separators for each sub-hc working from bottom of tree to top */
562:   /* this looks like nsep[]=segments */
563:   stages = (int *) bss_malloc((level+1)*INT_LEN);
564:   segs   = (int *) bss_malloc((level+1)*INT_LEN);
565:   ivec_zero(stages,level+1);
566:   ivec_copy(segs,nsep,level+1);
567:   for (i=0; i<level; i++)
568:     {segs[i+1] += segs[i];}
569:   stages[0] = segs[0];

571:   /* temporary vectors  */
572:   u  = (REAL *) bss_malloc(n*sizeof(REAL));
573:   z  = (REAL *) bss_malloc(n*sizeof(REAL));
574:   v  = (REAL *) bss_malloc(a_m*sizeof(REAL));
575:   uu = (REAL *) bss_malloc(m*sizeof(REAL));
576:   w  = (REAL *) bss_malloc(m*sizeof(REAL));

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

582:   /* storage for sparse x values */
583:   n_global = xxt_handle->info->n_global;
584:   xxt_max_nnz = (int)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/num_nodes;
585:   x = (REAL *) bss_malloc(xxt_max_nnz*sizeof(REAL));
586:   xxt_nnz = 0;

588:   /* LATER - can embed next sep to fire in gs */
589:   /* time to make the donuts - generate X factor */
590:   for (dim=i=j=0;i<m;i++)
591:     {
592:       /* time to move to the next level? */
593:       while (i==segs[dim])
594:         {
595: #ifdef SAFE          
596:           if (dim==level)
597:             {error_msg_fatal("dim about to exceed level\n"); break;}
598: #endif

600:           stages[dim++]=i;
601:           end+=lnsep[dim];
602:         }
603:       stages[dim]=i;

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

611:       /* shouldn't need this */
612:       if (col==INT_MAX)
613:         {
614:           error_msg_warning("hey ... col==INT_MAX??\n");
615:           continue;
616:         }

618:       /* do I own it? I should */
619:       rvec_zero(v ,a_m);
620:       if (col==fo[start])
621:         {
622:           start++;
623:           idex=ivec_linear_search(col, a_local2global, a_n);
624:           if (idex!=-1)
625:             {v[idex] = 1.0; j++;}
626:           else
627:             {error_msg_fatal("NOT FOUND!\n");}
628:         }
629:       else
630:         {
631:           idex=ivec_linear_search(col, a_local2global, a_m);
632:           if (idex!=-1)
633:             {v[idex] = 1.0;}
634:         }

636:       /* perform u = A.v_l */
637:       rvec_zero(u,n);
638:       do_matvec(xxt_handle->mvi,v,u);

640:       /* uu =  X^T.u_l (local portion) */
641:       /* technically only need to zero out first i entries */
642:       /* later turn this into an XXT_solve call ? */
643:       rvec_zero(uu,m);
644:       x_ptr=x;
645:       iptr = col_indices;
646:       for (k=0; k<i; k++)
647:         {
648:           off = *iptr++;
649:           len = *iptr++;

651: #if   BLAS||CBLAS
652:           uu[k] = dot(len,u+off,1,x_ptr,1);
653: #else
654:           uu[k] = rvec_dot(u+off,x_ptr,len);
655: #endif
656:           x_ptr+=len;
657:         }


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

663:       /* z = X.uu */
664:       rvec_zero(z,n);
665:       x_ptr=x;
666:       iptr = col_indices;
667:       for (k=0; k<i; k++)
668:         {
669:           off = *iptr++;
670:           len = *iptr++;

672: #if   BLAS||CBLAS
673:           axpy(len,uu[k],x_ptr,1,z+off,1);
674: #else
675:           rvec_axpy(z+off,x_ptr,uu[k],len);
676: #endif
677:           x_ptr+=len;
678:         }

680:       /* compute v_l = v_l - z */
681:       rvec_zero(v+a_n,a_m-a_n);
682: #if   BLAS&&CBLAS
683:       axpy(n,-1.0,z,1,v,1);
684: #else
685:       rvec_axpy(v,z,-1.0,n);
686: #endif

688:       /* compute u_l = A.v_l */
689:       if (a_n!=a_m)
690:         {gs_gop_hc(gs_handle,v,"+\0",dim);}
691:       rvec_zero(u,n);
692:       do_matvec(xxt_handle->mvi,v,u);

694:       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */
695: #if   BLAS&&CBLAS
696:       alpha = dot(n,u,1,v,1);
697: #else
698:       alpha = rvec_dot(u,v,n);
699: #endif
700:       /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */
701:       grop_hc(&alpha, &alpha_w, 1, op, dim);

703:       alpha = (REAL) sqrt((double)alpha);

705:       /* check for small alpha                             */
706:       /* LATER use this to detect and determine null space */
707: #ifdef r8
708:       if (fabs(alpha)<1.0e-14)
709:         {error_msg_fatal("bad alpha! %g\n",alpha);}
710: #else
711:       if (fabs((double) alpha) < 1.0e-6)
712:         {error_msg_fatal("bad alpha! %g\n",alpha);}
713: #endif

715:       /* compute v_l = v_l/sqrt(alpha) */
716:       rvec_scale(v,1.0/alpha,n);

718:       /* add newly generated column, v_l, to X */
719:       flag = 1;
720:       off=len=0;
721:       for (k=0; k<n; k++)
722:         {
723:           if (v[k]!=0.0)
724:             {
725:               len=k;
726:               if (flag)
727:                 {off=k; flag=0;}
728:             }
729:         }

731:       len -= (off-1);

733:       if (len>0)
734:         {
735:           if ((xxt_nnz+len)>xxt_max_nnz)
736:             {
737:               error_msg_warning("increasing space for X by 2x!\n");
738:               xxt_max_nnz *= 2;
739:               x_ptr = (REAL *) bss_malloc(xxt_max_nnz*sizeof(REAL));
740:               rvec_copy(x_ptr,x,xxt_nnz);
741:               bss_free(x);
742:               x = x_ptr;
743:               x_ptr+=xxt_nnz;
744:             }
745:           xxt_nnz += len;
746:           rvec_copy(x_ptr,v+off,len);

748:           /* keep track of number of zeros */
749:           if (dim)
750:             {
751:               for (k=0; k<len; k++)
752:                 {
753:                   if (x_ptr[k]==0.0)
754:                     {xxt_zero_nnz++;}
755:                 }
756:             }
757:           else
758:             {
759:               for (k=0; k<len; k++)
760:                 {
761:                   if (x_ptr[k]==0.0)
762:                     {xxt_zero_nnz_0++;}
763:                 }
764:             }
765:           col_indices[2*i] = off;
766:           col_sz[i] = col_indices[2*i+1] = len;
767:           col_vals[i] = x_ptr;
768:         }
769:       else
770:         {
771:           col_indices[2*i] = 0;
772:           col_sz[i] = col_indices[2*i+1] = 0;
773:           col_vals[i] = x_ptr;
774:         }
775:     }

777:   /* close off stages for execution phase */
778:   while (dim!=level)
779:     {
780:       stages[dim++]=i;
781:       error_msg_warning("disconnected!!! dim(%d)!=level(%d)\n",dim,level);
782:     }
783:   stages[dim]=i;

785:   xxt_handle->info->n=xxt_handle->mvi->n;
786:   xxt_handle->info->m=m;
787:   xxt_handle->info->nnz=xxt_nnz;
788:   xxt_handle->info->max_nnz=xxt_max_nnz;
789:   xxt_handle->info->msg_buf_sz=stages[level]-stages[0];
790:   xxt_handle->info->solve_uu = (REAL *) bss_malloc(m*sizeof(REAL));
791:   xxt_handle->info->solve_w  = (REAL *) bss_malloc(m*sizeof(REAL));
792:   xxt_handle->info->x=x;
793:   xxt_handle->info->col_vals=col_vals;
794:   xxt_handle->info->col_sz=col_sz;
795:   xxt_handle->info->col_indices=col_indices;
796:   xxt_handle->info->stages=stages;
797:   xxt_handle->info->nsolves=0;
798:   xxt_handle->info->tot_solve_time=0.0;

800:   bss_free(segs);
801:   bss_free(u);
802:   bss_free(v);
803:   bss_free(uu);
804:   bss_free(z);
805:   bss_free(w);

807: #ifdef DEBUG
808:   error_msg_warning("xxt_generate() :: end\n");
809: #endif

811:   return(0);
812: }


815: /*************************************xxt.c************************************
816: Function: 

818: Input : 
819: Output: 
820: Return: 
821: Description:  
822: **************************************xxt.c***********************************/
823: static
824: void
825: do_xxt_solve(xxt_ADT xxt_handle, register REAL *uc)
826: {
827:   register int off, len, *iptr;
828:   int level       =xxt_handle->level;
829:   int n           =xxt_handle->info->n;
830:   int m           =xxt_handle->info->m;
831:   int *stages     =xxt_handle->info->stages;
832:   int *col_indices=xxt_handle->info->col_indices;
833:   register REAL *x_ptr, *uu_ptr;
834: #if   BLAS||CBLAS
835:   REAL zero=0.0;
836: #endif
837:   REAL *solve_uu=xxt_handle->info->solve_uu;
838:   REAL *solve_w =xxt_handle->info->solve_w;
839:   REAL *x       =xxt_handle->info->x;

841: #ifdef DEBUG
842:   error_msg_warning("do_xxt_solve() :: begin\n");
843: #endif

845:   uu_ptr=solve_uu;
846: #if   BLAS||CBLAS
847:   copy(m,&zero,0,uu_ptr,1);
848: #else
849:   rvec_zero(uu_ptr,m);
850: #endif

852:   /* x  = X.Y^T.b */
853:   /* uu = Y^T.b */
854:   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
855:     {
856:       off=*iptr++; len=*iptr++;
857: #if   BLAS||CBLAS
858:       *uu_ptr++ = dot(len,uc+off,1,x_ptr,1);
859: #else
860:       *uu_ptr++ = rvec_dot(uc+off,x_ptr,len);
861: #endif
862:     }

864:   /* comunication of beta */
865:   uu_ptr=solve_uu;
866:   if (level) {ssgl_radd(uu_ptr, solve_w, level, stages);}

868: #if   BLAS||CBLAS
869:   copy(n,&zero,0,uc,1);
870: #else
871:   rvec_zero(uc,n);
872: #endif

874:   /* x = X.uu */
875:   for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len)
876:     {
877:       off=*iptr++; len=*iptr++;
878: #if   BLAS||CBLAS
879:       axpy(len,*uu_ptr++,x_ptr,1,uc+off,1);
880: #else
881:       rvec_axpy(uc+off,x_ptr,*uu_ptr++,len);
882: #endif
883:     }

885: #ifdef DEBUG
886:   error_msg_warning("do_xxt_solve() :: end\n");
887: #endif
888: }


891: /*************************************Xxt.c************************************
892: Function: check_init

894: Input :
895: Output:
896: Return:
897: Description:
898: **************************************xxt.c***********************************/
899: static
900: void
901: check_init(void)
902: {
903: #ifdef DEBUG
904:   error_msg_warning("check_init() :: start %d\n",n_xxt_handles);
905: #endif

907:   comm_init();
908:   /*
909:   perm_init(); 
910:   bss_init();
911:   */

913: #ifdef DEBUG
914:   error_msg_warning("check_init() :: end   %d\n",n_xxt_handles);
915: #endif
916: }


919: /*************************************xxt.c************************************
920: Function: check_handle()

922: Input :
923: Output:
924: Return:
925: Description:
926: **************************************xxt.c***********************************/
927: static
928: void 
929: check_handle(xxt_ADT xxt_handle)
930: {
931: #ifdef SAFE
932:   int vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
933: #endif


936: #ifdef DEBUG
937:   error_msg_warning("check_handle() :: start %d\n",n_xxt_handles);
938: #endif

940:   if (xxt_handle==NULL)
941:     {error_msg_fatal("check_handle() :: bad handle :: NULL %d\n",xxt_handle);}

943: #ifdef SAFE
944:   vals[0]=vals[1]=xxt_handle->id;
945:   giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
946:   if ((vals[0]!=vals[1])||(xxt_handle->id<=0))
947:     {error_msg_fatal("check_handle() :: bad handle :: id mismatch min/max %d/%d %d\n",
948:                      vals[0],vals[1], xxt_handle->id);}
949: #endif

951: #ifdef DEBUG
952:   error_msg_warning("check_handle() :: end   %d\n",n_xxt_handles);
953: #endif
954: }


957: /*************************************xxt.c************************************
958: Function: det_separators

960: Input :
961: Output:
962: Return:
963: Description:
964:   det_separators(xxt_handle, local2global, n, m, mylocmatvec, grid_data);
965: **************************************xxt.c***********************************/
966: static 
967: void 
968: det_separators(xxt_ADT xxt_handle)
969: {
970:   int i, ct, id;
971:   int mask, edge, *iptr;
972:   int *dir, *used;
973:   int sum[4], w[4];
974:   REAL rsum[4], rw[4];
975:   int op[] = {GL_ADD,0};
976:   REAL *lhs, *rhs;
977:   int *nsep, *lnsep, *fo, nfo=0;
978:   gs_ADT gs_handle=xxt_handle->mvi->gs_handle;
979:   int *local2global=xxt_handle->mvi->local2global;
980:   int  n=xxt_handle->mvi->n;
981:   int  m=xxt_handle->mvi->m;
982:   int level=xxt_handle->level;
983:   int shared=FALSE;

985: #ifdef DEBUG
986:   error_msg_warning("det_separators() :: start %d %d %d\n",level,n,m);
987: #endif
988: 
989:   dir  = (int*)bss_malloc(INT_LEN*(level+1));
990:   nsep = (int*)bss_malloc(INT_LEN*(level+1));
991:   lnsep= (int*)bss_malloc(INT_LEN*(level+1));
992:   fo   = (int*)bss_malloc(INT_LEN*(n+1));
993:   used = (int*)bss_malloc(INT_LEN*n);

995:   ivec_zero(dir  ,level+1);
996:   ivec_zero(nsep ,level+1);
997:   ivec_zero(lnsep,level+1);
998:   ivec_set (fo   ,-1,n+1);
999:   ivec_zero(used,n);

1001:   lhs  = (double*)bss_malloc(REAL_LEN*m);
1002:   rhs  = (double*)bss_malloc(REAL_LEN*m);

1004:   /* determine the # of unique dof */
1005:   rvec_zero(lhs,m);
1006:   rvec_set(lhs,1.0,n);
1007:   gs_gop_hc(gs_handle,lhs,"+\0",level);
1008:   rvec_zero(rsum,2);
1009:   for (ct=i=0;i<n;i++)
1010:     {
1011:       if (lhs[i]!=0.0)
1012:         {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];}
1013:     }
1014:   grop_hc(rsum,rw,2,op,level);
1015:   rsum[0]+=0.1;
1016:   rsum[1]+=0.1;
1017:   /*  if (!my_id)
1018:     {
1019:       printf("xxt n unique = %d (%g)\n",(int) rsum[0], rsum[0]);
1020:       printf("xxt n shared = %d (%g)\n",(int) rsum[1], rsum[1]);
1021:       }*/

1023:   if (fabs(rsum[0]-rsum[1])>EPS)
1024:     {shared=TRUE;}

1026:   xxt_handle->info->n_global=xxt_handle->info->m_global=(int) rsum[0];
1027:   xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(int) rsum[0];

1029:   /* determine separator sets top down */
1030:   if (shared)
1031:     {
1032:       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
1033:         {
1034:           /* set rsh of hc, fire, and collect lhs responses */
1035:           (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
1036:           gs_gop_hc(gs_handle,lhs,"+\0",edge);
1037: 
1038:           /* set lsh of hc, fire, and collect rhs responses */
1039:           (id<mask) ? rvec_set(rhs,1.0,m) : rvec_zero(rhs,m);
1040:           gs_gop_hc(gs_handle,rhs,"+\0",edge);
1041: 
1042:           for (i=0;i<n;i++)
1043:             {
1044:               if (id< mask)
1045:                 {
1046:                   if (lhs[i]!=0.0)
1047:                     {lhs[i]=1.0;}
1048:                 }
1049:               if (id>=mask)
1050:                 {
1051:                   if (rhs[i]!=0.0)
1052:                     {rhs[i]=1.0;}
1053:                 }
1054:             }

1056:           if (id< mask)
1057:             {gs_gop_hc(gs_handle,lhs,"+\0",edge-1);}
1058:           else
1059:             {gs_gop_hc(gs_handle,rhs,"+\0",edge-1);}

1061:           /* count number of dofs I own that have signal and not in sep set */
1062:           rvec_zero(rsum,4);
1063:           for (ivec_zero(sum,4),ct=i=0;i<n;i++)
1064:             {
1065:               if (!used[i])
1066:                 {
1067:                   /* number of unmarked dofs on node */
1068:                   ct++;
1069:                   /* number of dofs to be marked on lhs hc */
1070:                   if (id< mask)
1071:                     {
1072:                       if (lhs[i]!=0.0)
1073:                         {sum[0]++; rsum[0]+=1.0/lhs[i];}
1074:                     }
1075:                   /* number of dofs to be marked on rhs hc */
1076:                   if (id>=mask)
1077:                     {
1078:                       if (rhs[i]!=0.0)
1079:                         {sum[1]++; rsum[1]+=1.0/rhs[i];}
1080:                     }
1081:                 }
1082:             }

1084:           /* go for load balance - choose half with most unmarked dofs, bias LHS */
1085:           (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
1086:           (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct);
1087:           giop_hc(sum,w,4,op,edge);
1088:           grop_hc(rsum,rw,4,op,edge);
1089:           rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1;

1091:           if (id<mask)
1092:             {
1093:               /* mark dofs I own that have signal and not in sep set */
1094:               for (ct=i=0;i<n;i++)
1095:                 {
1096:                   if ((!used[i])&&(lhs[i]!=0.0))
1097:                     {
1098:                       ct++; nfo++;

1100:                       if (nfo>n)
1101:                         {error_msg_fatal("nfo about to exceed n\n");}

1103:                       *--iptr = local2global[i];
1104:                       used[i]=edge;
1105:                     }
1106:                 }
1107:               if (ct>1) {ivec_sort(iptr,ct);}

1109:               lnsep[edge]=ct;
1110:               nsep[edge]=(int) rsum[0];
1111:               dir [edge]=LEFT;
1112:             }

1114:           if (id>=mask)
1115:             {
1116:               /* mark dofs I own that have signal and not in sep set */
1117:               for (ct=i=0;i<n;i++)
1118:                 {
1119:                   if ((!used[i])&&(rhs[i]!=0.0))
1120:                     {
1121:                       ct++; nfo++;

1123:                       if (nfo>n)
1124:                         {error_msg_fatal("nfo about to exceed n\n");}

1126:                       *--iptr = local2global[i];
1127:                       used[i]=edge;
1128:                     }
1129:                 }
1130:               if (ct>1) {ivec_sort(iptr,ct);}

1132:               lnsep[edge]=ct;
1133:               nsep[edge]= (int) rsum[1];
1134:               dir [edge]=RIGHT;
1135:             }

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

1140:           /* fold rhs hc into lower */
1141:           if (id>=mask)
1142:             {id-=mask;}
1143:         }
1144:     }
1145:   else
1146:     {
1147:       for (iptr=fo+n,id=my_id,mask=num_nodes>>1,edge=level;edge>0;edge--,mask>>=1)
1148:         {
1149:           /* set rsh of hc, fire, and collect lhs responses */
1150:           (id<mask) ? rvec_zero(lhs,m) : rvec_set(lhs,1.0,m);
1151:           gs_gop_hc(gs_handle,lhs,"+\0",edge);

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

1157:           /* count number of dofs I own that have signal and not in sep set */
1158:           for (ivec_zero(sum,4),ct=i=0;i<n;i++)
1159:             {
1160:               if (!used[i])
1161:                 {
1162:                   /* number of unmarked dofs on node */
1163:                   ct++;
1164:                   /* number of dofs to be marked on lhs hc */
1165:                   if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;}
1166:                   /* number of dofs to be marked on rhs hc */
1167:                   if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;}
1168:                 }
1169:             }

1171:           /* go for load balance - choose half with most unmarked dofs, bias LHS */
1172:           (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
1173:           giop_hc(sum,w,4,op,edge);

1175:           /* lhs hc wins */
1176:           if (sum[2]>=sum[3])
1177:             {
1178:               if (id<mask)
1179:                 {
1180:                   /* mark dofs I own that have signal and not in sep set */
1181:                   for (ct=i=0;i<n;i++)
1182:                     {
1183:                       if ((!used[i])&&(lhs[i]!=0.0))
1184:                         {
1185:                           ct++; nfo++;
1186:                           *--iptr = local2global[i];
1187:                           used[i]=edge;
1188:                         }
1189:                     }
1190:                   if (ct>1) {ivec_sort(iptr,ct);}
1191:                   lnsep[edge]=ct;
1192:                 }
1193:               nsep[edge]=sum[0];
1194:               dir [edge]=LEFT;
1195:             }
1196:           /* rhs hc wins */
1197:           else
1198:             {
1199:               if (id>=mask)
1200:                 {
1201:                   /* mark dofs I own that have signal and not in sep set */
1202:                   for (ct=i=0;i<n;i++)
1203:                     {
1204:                       if ((!used[i])&&(rhs[i]!=0.0))
1205:                         {
1206:                           ct++; nfo++;
1207:                           *--iptr = local2global[i];
1208:                           used[i]=edge;
1209:                         }
1210:                     }
1211:                   if (ct>1) {ivec_sort(iptr,ct);}
1212:                   lnsep[edge]=ct;
1213:                 }
1214:               nsep[edge]=sum[1];
1215:               dir [edge]=RIGHT;
1216:             }
1217:           /* LATER or we can recur on these to order seps at this level */
1218:           /* do we need full set of separators for this?                */

1220:           /* fold rhs hc into lower */
1221:           if (id>=mask)
1222:             {id-=mask;}
1223:         }
1224:     }


1227:   /* level 0 is on processor case - so mark the remainder */
1228:   for (ct=i=0;i<n;i++)
1229:     {
1230:       if (!used[i])
1231:         {
1232:           ct++; nfo++;
1233:           *--iptr = local2global[i];
1234:           used[i]=edge;
1235:         }
1236:     }
1237:   if (ct>1) {ivec_sort(iptr,ct);}
1238:   lnsep[edge]=ct;
1239:   nsep [edge]=ct;
1240:   dir  [edge]=LEFT;

1242:   xxt_handle->info->nsep=nsep;
1243:   xxt_handle->info->lnsep=lnsep;
1244:   xxt_handle->info->fo=fo;
1245:   xxt_handle->info->nfo=nfo;

1247:   bss_free(dir);
1248:   bss_free(lhs);
1249:   bss_free(rhs);
1250:   bss_free(used);

1252: #ifdef DEBUG  
1253:   error_msg_warning("det_separators() :: end\n");
1254: #endif
1255: }


1258: /*************************************xxt.c************************************
1259: Function: set_mvi

1261: Input :
1262: Output:
1263: Return:
1264: Description:
1265: **************************************xxt.c***********************************/
1266: static
1267: mv_info *set_mvi(int *local2global, int n, int m, void *matvec, void *grid_data)
1268: {
1269:   mv_info *mvi;


1272: #ifdef DEBUG
1273:   error_msg_warning("set_mvi() :: start\n");
1274: #endif

1276:   mvi = (mv_info*)bss_malloc(sizeof(mv_info));
1277:   mvi->n=n;
1278:   mvi->m=m;
1279:   mvi->n_global=-1;
1280:   mvi->m_global=-1;
1281:   mvi->local2global=(int*)bss_malloc((m+1)*INT_LEN);
1282:   ivec_copy(mvi->local2global,local2global,m);
1283:   mvi->local2global[m] = INT_MAX;
1284:   mvi->matvec=(int (*)(mv_info*,REAL*,REAL*))matvec;
1285:   mvi->grid_data=grid_data;

1287:   /* set xxt communication handle to perform restricted matvec */
1288:   mvi->gs_handle = gs_init(local2global, m, num_nodes);

1290: #ifdef DEBUG
1291:   error_msg_warning("set_mvi() :: end   \n");
1292: #endif
1293: 
1294:   return(mvi);
1295: }


1298: /*************************************xxt.c************************************
1299: Function: set_mvi

1301: Input :
1302: Output:
1303: Return:
1304: Description:

1306:       computes u = A.v 
1307:       do_matvec(xxt_handle->mvi,v,u);
1308: **************************************xxt.c***********************************/
1309: static
1310: void do_matvec(mv_info *A, REAL *v, REAL *u)
1311: {
1312:   A->matvec((mv_info*)A->grid_data,v,u);
1313: }