Actual source code: xyt.c
2: /*************************************xyt.c************************************
3: Module Name: xyt
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: **************************************xyt.c***********************************/
19: #include <../src/ksp/pc/impls/tfs/tfs.h>
21: #define LEFT -1
22: #define RIGHT 1
23: #define BOTH 0
25: typedef struct xyt_solver_info {
26: PetscInt n, m, n_global, m_global;
27: PetscInt nnz, max_nnz, msg_buf_sz;
28: PetscInt *nsep, *lnsep, *fo, nfo, *stages;
29: PetscInt *xcol_sz, *xcol_indices;
30: PetscScalar **xcol_vals, *x, *solve_uu, *solve_w;
31: PetscInt *ycol_sz, *ycol_indices;
32: PetscScalar **ycol_vals, *y;
33: PetscInt nsolves;
34: PetscScalar tot_solve_time;
35: } xyt_info;
37: typedef struct matvec_info {
38: PetscInt n, m, n_global, m_global;
39: PetscInt *local2global;
40: PCTFS_gs_ADT PCTFS_gs_handle;
41: PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*);
42: void *grid_data;
43: } mv_info;
45: struct xyt_CDT {
46: PetscInt id;
47: PetscInt ns;
48: PetscInt level;
49: xyt_info *info;
50: mv_info *mvi;
51: };
53: static PetscInt n_xyt =0;
54: static PetscInt n_xyt_handles=0;
56: /* prototypes */
57: static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs);
58: static PetscErrorCode check_handle(xyt_ADT xyt_handle);
59: static PetscErrorCode det_separators(xyt_ADT xyt_handle);
60: static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u);
61: static PetscErrorCode xyt_generate(xyt_ADT xyt_handle);
62: static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle);
63: static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data);
65: /**************************************xyt.c***********************************/
66: xyt_ADT XYT_new(void)
67: {
68: xyt_ADT xyt_handle;
70: /* rolling count on n_xyt ... pot. problem here */
71: n_xyt_handles++;
72: xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT));
73: xyt_handle->id = ++n_xyt;
74: xyt_handle->info = NULL;
75: xyt_handle->mvi = NULL;
77: return(xyt_handle);
78: }
80: /**************************************xyt.c***********************************/
81: PetscErrorCode XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */
82: PetscInt *local2global, /* global column mapping */
83: PetscInt n, /* local num rows */
84: PetscInt m, /* local num cols */
85: PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */
86: void *grid_data) /* grid data for matvec */
87: {
89: PCTFS_comm_init();
90: check_handle(xyt_handle);
92: /* only 2^k for now and all nodes participating */
95: /* space for X info */
96: xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info));
98: /* set up matvec handles */
99: xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data);
101: /* matrix is assumed to be of full rank */
102: /* LATER we can reset to indicate rank def. */
103: xyt_handle->ns=0;
105: /* determine separators and generate firing order - NB xyt info set here */
106: det_separators(xyt_handle);
108: return(do_xyt_factor(xyt_handle));
109: }
111: /**************************************xyt.c***********************************/
112: PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b)
113: {
114: PCTFS_comm_init();
115: check_handle(xyt_handle);
117: /* need to copy b into x? */
118: if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n);
119: return do_xyt_solve(xyt_handle,x);
120: }
122: /**************************************xyt.c***********************************/
123: PetscErrorCode XYT_free(xyt_ADT xyt_handle)
124: {
125: PCTFS_comm_init();
126: check_handle(xyt_handle);
127: n_xyt_handles--;
129: free(xyt_handle->info->nsep);
130: free(xyt_handle->info->lnsep);
131: free(xyt_handle->info->fo);
132: free(xyt_handle->info->stages);
133: free(xyt_handle->info->solve_uu);
134: free(xyt_handle->info->solve_w);
135: free(xyt_handle->info->x);
136: free(xyt_handle->info->xcol_vals);
137: free(xyt_handle->info->xcol_sz);
138: free(xyt_handle->info->xcol_indices);
139: free(xyt_handle->info->y);
140: free(xyt_handle->info->ycol_vals);
141: free(xyt_handle->info->ycol_sz);
142: free(xyt_handle->info->ycol_indices);
143: free(xyt_handle->info);
144: free(xyt_handle->mvi->local2global);
145: PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle);
146: free(xyt_handle->mvi);
147: free(xyt_handle);
149: /* if the check fails we nuke */
150: /* if NULL pointer passed to free we nuke */
151: /* if the calls to free fail that's not my problem */
152: return(0);
153: }
155: /**************************************xyt.c***********************************/
156: PetscErrorCode XYT_stats(xyt_ADT xyt_handle)
157: {
158: PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD};
159: PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD};
160: PetscInt vals[9], work[9];
161: PetscScalar fvals[3], fwork[3];
163: PCTFS_comm_init();
164: check_handle(xyt_handle);
166: /* if factorization not done there are no stats */
167: if (!xyt_handle->info||!xyt_handle->mvi) {
168: if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n");
169: return 1;
170: }
172: vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz;
173: vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n;
174: vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz;
175: PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
177: fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++;
178: PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop);
180: if (!PCTFS_my_id) {
181: PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_nnz=%D\n",PCTFS_my_id,vals[0]);
182: PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_nnz=%D\n",PCTFS_my_id,vals[1]);
183: PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);
184: PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_nnz=%D\n",PCTFS_my_id,vals[2]);
185: PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(2d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5)));
186: PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(3d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667)));
187: PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_n =%D\n",PCTFS_my_id,vals[3]);
188: PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_n =%D\n",PCTFS_my_id,vals[4]);
189: PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);
190: PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_n =%D\n",PCTFS_my_id,vals[5]);
191: PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_buf=%D\n",PCTFS_my_id,vals[6]);
192: PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_buf=%D\n",PCTFS_my_id,vals[7]);
193: PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);
194: PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_slv=%g\n",PCTFS_my_id,fvals[0]);
195: PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_slv=%g\n",PCTFS_my_id,fvals[1]);
196: PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);
197: }
199: return(0);
200: }
202: /*************************************xyt.c************************************
204: Description: get A_local, local portion of global coarse matrix which
205: is a row dist. nxm matrix w/ n<m.
206: o my_ml holds address of ML struct associated w/A_local and coarse grid
207: o local2global holds global number of column i (i=0,...,m-1)
208: o local2global holds global number of row i (i=0,...,n-1)
209: o mylocmatvec performs A_local . vec_local (note that gs is performed using
210: PCTFS_gs_init/gop).
212: mylocmatvec = my_ml->Amat[grid_tag].matvec->external;
213: mylocmatvec (void :: void *data, double *in, double *out)
214: **************************************xyt.c***********************************/
215: static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle)
216: {
217: return xyt_generate(xyt_handle);
218: }
220: /**************************************xyt.c***********************************/
221: static PetscErrorCode xyt_generate(xyt_ADT xyt_handle)
222: {
223: PetscInt i,j,k,idx;
224: PetscInt dim, col;
225: PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w;
226: PetscInt *segs;
227: PetscInt op[] = {GL_ADD,0};
228: PetscInt off, len;
229: PetscScalar *x_ptr, *y_ptr;
230: PetscInt *iptr, flag;
231: PetscInt start =0, end, work;
232: PetscInt op2[] = {GL_MIN,0};
233: PCTFS_gs_ADT PCTFS_gs_handle;
234: PetscInt *nsep, *lnsep, *fo;
235: PetscInt a_n =xyt_handle->mvi->n;
236: PetscInt a_m =xyt_handle->mvi->m;
237: PetscInt *a_local2global=xyt_handle->mvi->local2global;
238: PetscInt level;
239: PetscInt n, m;
240: PetscInt *xcol_sz, *xcol_indices, *stages;
241: PetscScalar **xcol_vals, *x;
242: PetscInt *ycol_sz, *ycol_indices;
243: PetscScalar **ycol_vals, *y;
244: PetscInt n_global;
245: PetscInt xt_nnz =0, xt_max_nnz=0;
246: PetscInt yt_nnz =0, yt_max_nnz=0;
247: PetscInt xt_zero_nnz =0;
248: PetscInt xt_zero_nnz_0=0;
249: PetscInt yt_zero_nnz =0;
250: PetscInt yt_zero_nnz_0=0;
251: PetscBLASInt i1 = 1,dlen;
252: PetscScalar dm1 = -1.0;
254: n =xyt_handle->mvi->n;
255: nsep =xyt_handle->info->nsep;
256: lnsep =xyt_handle->info->lnsep;
257: fo =xyt_handle->info->fo;
258: end =lnsep[0];
259: level =xyt_handle->level;
260: PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
262: /* is there a null space? */
263: /* LATER add in ability to detect null space by checking alpha */
264: for (i=0, j=0; i<=level; i++) j+=nsep[i];
266: m = j-xyt_handle->ns;
267: if (m!=j) {
268: PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);
269: }
271: PetscInfo(0,"xyt_generate() :: X(%D,%D)\n",n,m);
273: /* get and initialize storage for x local */
274: /* note that x local is nxm and stored by columns */
275: xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
276: xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
277: xcol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
278: for (i=j=0; i<m; i++, j+=2) {
279: xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1;
280: xcol_vals[i] = NULL;
281: }
282: xcol_indices[j]=-1;
284: /* get and initialize storage for y local */
285: /* note that y local is nxm and stored by columns */
286: ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt));
287: ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt));
288: ycol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*));
289: for (i=j=0; i<m; i++, j+=2) {
290: ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1;
291: ycol_vals[i] = NULL;
292: }
293: ycol_indices[j]=-1;
295: /* size of separators for each sub-hc working from bottom of tree to top */
296: /* this looks like nsep[]=segments */
297: stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
298: segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt));
299: PCTFS_ivec_zero(stages,level+1);
300: PCTFS_ivec_copy(segs,nsep,level+1);
301: for (i=0; i<level; i++) segs[i+1] += segs[i];
302: stages[0] = segs[0];
304: /* temporary vectors */
305: u = (PetscScalar*) malloc(n*sizeof(PetscScalar));
306: z = (PetscScalar*) malloc(n*sizeof(PetscScalar));
307: v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar));
308: uu = (PetscScalar*) malloc(m*sizeof(PetscScalar));
309: w = (PetscScalar*) malloc(m*sizeof(PetscScalar));
311: /* extra nnz due to replication of vertices across separators */
312: for (i=1, j=0; i<=level; i++) j+=nsep[i];
314: /* storage for sparse x values */
315: n_global = xyt_handle->info->n_global;
316: xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes;
317: x = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar));
318: y = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar));
320: /* LATER - can embed next sep to fire in gs */
321: /* time to make the donuts - generate X factor */
322: for (dim=i=j=0; i<m; i++) {
323: /* time to move to the next level? */
324: while (i==segs[dim]) {
326: stages[dim++]=i;
327: end +=lnsep[dim];
328: }
329: stages[dim]=i;
331: /* which column are we firing? */
332: /* i.e. set v_l */
333: /* use new seps and do global min across hc to determine which one to fire */
334: (start<end) ? (col=fo[start]) : (col=INT_MAX);
335: PCTFS_giop_hc(&col,&work,1,op2,dim);
337: /* shouldn't need this */
338: if (col==INT_MAX) {
339: PetscInfo(0,"hey ... col==INT_MAX??\n");
340: continue;
341: }
343: /* do I own it? I should */
344: PCTFS_rvec_zero(v,a_m);
345: if (col==fo[start]) {
346: start++;
347: idx=PCTFS_ivec_linear_search(col, a_local2global, a_n);
348: if (idx!=-1) {
349: v[idx] = 1.0;
350: j++;
351: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!");
352: } else {
353: idx=PCTFS_ivec_linear_search(col, a_local2global, a_m);
354: if (idx!=-1) v[idx] = 1.0;
355: }
357: /* perform u = A.v_l */
358: PCTFS_rvec_zero(u,n);
359: do_matvec(xyt_handle->mvi,v,u);
361: /* uu = X^T.u_l (local portion) */
362: /* technically only need to zero out first i entries */
363: /* later turn this into an XYT_solve call ? */
364: PCTFS_rvec_zero(uu,m);
365: y_ptr=y;
366: iptr = ycol_indices;
367: for (k=0; k<i; k++) {
368: off = *iptr++;
369: len = *iptr++;
370: PetscBLASIntCast(len,&dlen);
371: PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1));
372: y_ptr+=len;
373: }
375: /* uu = X^T.u_l (comm portion) */
376: PCTFS_ssgl_radd (uu, w, dim, stages);
378: /* z = X.uu */
379: PCTFS_rvec_zero(z,n);
380: x_ptr=x;
381: iptr = xcol_indices;
382: for (k=0; k<i; k++) {
383: off = *iptr++;
384: len = *iptr++;
385: PetscBLASIntCast(len,&dlen);
386: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1));
387: x_ptr+=len;
388: }
390: /* compute v_l = v_l - z */
391: PCTFS_rvec_zero(v+a_n,a_m-a_n);
392: PetscBLASIntCast(n,&dlen);
393: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1));
395: /* compute u_l = A.v_l */
396: if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim);
397: PCTFS_rvec_zero(u,n);
398: do_matvec(xyt_handle->mvi,v,u);
400: /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */
401: PetscBLASIntCast(n,&dlen);
402: PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,u,&i1));
403: /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */
404: PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim);
406: alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha);
408: /* check for small alpha */
409: /* LATER use this to detect and determine null space */
412: /* compute v_l = v_l/sqrt(alpha) */
413: PCTFS_rvec_scale(v,1.0/alpha,n);
414: PCTFS_rvec_scale(u,1.0/alpha,n);
416: /* add newly generated column, v_l, to X */
417: flag = 1;
418: off =len=0;
419: for (k=0; k<n; k++) {
420: if (v[k]!=0.0) {
421: len=k;
422: if (flag) {off=k; flag=0;}
423: }
424: }
426: len -= (off-1);
428: if (len>0) {
429: if ((xt_nnz+len)>xt_max_nnz) {
430: PetscInfo(0,"increasing space for X by 2x!\n");
431: xt_max_nnz *= 2;
432: x_ptr = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar));
433: PCTFS_rvec_copy(x_ptr,x,xt_nnz);
434: free(x);
435: x = x_ptr;
436: x_ptr+=xt_nnz;
437: }
438: xt_nnz += len;
439: PCTFS_rvec_copy(x_ptr,v+off,len);
441: /* keep track of number of zeros */
442: if (dim) {
443: for (k=0; k<len; k++) {
444: if (x_ptr[k]==0.0) xt_zero_nnz++;
445: }
446: } else {
447: for (k=0; k<len; k++) {
448: if (x_ptr[k]==0.0) xt_zero_nnz_0++;
449: }
450: }
451: xcol_indices[2*i] = off;
452: xcol_sz[i] = xcol_indices[2*i+1] = len;
453: xcol_vals[i] = x_ptr;
454: } else {
455: xcol_indices[2*i] = 0;
456: xcol_sz[i] = xcol_indices[2*i+1] = 0;
457: xcol_vals[i] = x_ptr;
458: }
460: /* add newly generated column, u_l, to Y */
461: flag = 1;
462: off =len=0;
463: for (k=0; k<n; k++) {
464: if (u[k]!=0.0) {
465: len=k;
466: if (flag) { off=k; flag=0; }
467: }
468: }
470: len -= (off-1);
472: if (len>0) {
473: if ((yt_nnz+len)>yt_max_nnz) {
474: PetscInfo(0,"increasing space for Y by 2x!\n");
475: yt_max_nnz *= 2;
476: y_ptr = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar));
477: PCTFS_rvec_copy(y_ptr,y,yt_nnz);
478: free(y);
479: y = y_ptr;
480: y_ptr+=yt_nnz;
481: }
482: yt_nnz += len;
483: PCTFS_rvec_copy(y_ptr,u+off,len);
485: /* keep track of number of zeros */
486: if (dim) {
487: for (k=0; k<len; k++) {
488: if (y_ptr[k]==0.0) yt_zero_nnz++;
489: }
490: } else {
491: for (k=0; k<len; k++) {
492: if (y_ptr[k]==0.0) yt_zero_nnz_0++;
493: }
494: }
495: ycol_indices[2*i] = off;
496: ycol_sz[i] = ycol_indices[2*i+1] = len;
497: ycol_vals[i] = y_ptr;
498: } else {
499: ycol_indices[2*i] = 0;
500: ycol_sz[i] = ycol_indices[2*i+1] = 0;
501: ycol_vals[i] = y_ptr;
502: }
503: }
505: /* close off stages for execution phase */
506: while (dim!=level) {
507: stages[dim++]=i;
508: PetscInfo(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);
509: }
510: stages[dim]=i;
512: xyt_handle->info->n =xyt_handle->mvi->n;
513: xyt_handle->info->m =m;
514: xyt_handle->info->nnz =xt_nnz + yt_nnz;
515: xyt_handle->info->max_nnz =xt_max_nnz + yt_max_nnz;
516: xyt_handle->info->msg_buf_sz =stages[level]-stages[0];
517: xyt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar));
518: xyt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar));
519: xyt_handle->info->x =x;
520: xyt_handle->info->xcol_vals =xcol_vals;
521: xyt_handle->info->xcol_sz =xcol_sz;
522: xyt_handle->info->xcol_indices=xcol_indices;
523: xyt_handle->info->stages =stages;
524: xyt_handle->info->y =y;
525: xyt_handle->info->ycol_vals =ycol_vals;
526: xyt_handle->info->ycol_sz =ycol_sz;
527: xyt_handle->info->ycol_indices=ycol_indices;
529: free(segs);
530: free(u);
531: free(v);
532: free(uu);
533: free(z);
534: free(w);
536: return(0);
537: }
539: /**************************************xyt.c***********************************/
540: static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc)
541: {
542: PetscInt off, len, *iptr;
543: PetscInt level =xyt_handle->level;
544: PetscInt n =xyt_handle->info->n;
545: PetscInt m =xyt_handle->info->m;
546: PetscInt *stages =xyt_handle->info->stages;
547: PetscInt *xcol_indices=xyt_handle->info->xcol_indices;
548: PetscInt *ycol_indices=xyt_handle->info->ycol_indices;
549: PetscScalar *x_ptr, *y_ptr, *uu_ptr;
550: PetscScalar *solve_uu=xyt_handle->info->solve_uu;
551: PetscScalar *solve_w =xyt_handle->info->solve_w;
552: PetscScalar *x =xyt_handle->info->x;
553: PetscScalar *y =xyt_handle->info->y;
554: PetscBLASInt i1 = 1,dlen;
556: uu_ptr=solve_uu;
557: PCTFS_rvec_zero(uu_ptr,m);
559: /* x = X.Y^T.b */
560: /* uu = Y^T.b */
561: for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) {
562: off =*iptr++;
563: len =*iptr++;
564: PetscBLASIntCast(len,&dlen);
565: PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1));
566: }
568: /* comunication of beta */
569: uu_ptr=solve_uu;
570: if (level) PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);
571: PCTFS_rvec_zero(uc,n);
573: /* x = X.uu */
574: for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) {
575: off =*iptr++;
576: len =*iptr++;
577: PetscBLASIntCast(len,&dlen);
578: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1));
579: }
580: return 0;
581: }
583: /**************************************xyt.c***********************************/
584: static PetscErrorCode check_handle(xyt_ADT xyt_handle)
585: {
586: PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX};
590: vals[0]=vals[1]=xyt_handle->id;
591: PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op);
593: return 0;
594: }
596: /**************************************xyt.c***********************************/
597: static PetscErrorCode det_separators(xyt_ADT xyt_handle)
598: {
599: PetscInt i, ct, id;
600: PetscInt mask, edge, *iptr;
601: PetscInt *dir, *used;
602: PetscInt sum[4], w[4];
603: PetscScalar rsum[4], rw[4];
604: PetscInt op[] = {GL_ADD,0};
605: PetscScalar *lhs, *rhs;
606: PetscInt *nsep, *lnsep, *fo, nfo=0;
607: PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle;
608: PetscInt *local2global =xyt_handle->mvi->local2global;
609: PetscInt n =xyt_handle->mvi->n;
610: PetscInt m =xyt_handle->mvi->m;
611: PetscInt level =xyt_handle->level;
612: PetscInt shared =0;
614: dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
615: nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
616: lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1));
617: fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1));
618: used = (PetscInt*)malloc(sizeof(PetscInt)*n);
620: PCTFS_ivec_zero(dir,level+1);
621: PCTFS_ivec_zero(nsep,level+1);
622: PCTFS_ivec_zero(lnsep,level+1);
623: PCTFS_ivec_set (fo,-1,n+1);
624: PCTFS_ivec_zero(used,n);
626: lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
627: rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m);
629: /* determine the # of unique dof */
630: PCTFS_rvec_zero(lhs,m);
631: PCTFS_rvec_set(lhs,1.0,n);
632: PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level);
633: PetscInfo(0,"done first PCTFS_gs_gop_hc\n");
634: PCTFS_rvec_zero(rsum,2);
635: for (i=0; i<n; i++) {
636: if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; }
637: if (lhs[i]!=1.0) shared=1;
638: }
640: PCTFS_grop_hc(rsum,rw,2,op,level);
641: rsum[0]+=0.1;
642: rsum[1]+=0.1;
644: xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0];
645: xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0];
647: /* determine separator sets top down */
648: if (shared) {
649: /* solution is to do as in the symmetric shared case but then */
650: /* pick the sub-hc with the most free dofs and do a mat-vec */
651: /* and pick up the responses on the other sub-hc from the */
652: /* initial separator set obtained from the symm. shared case */
653: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!");
654: /* [dead code deleted since it is unlikely to be completed] */
655: } else {
656: for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) {
657: /* set rsh of hc, fire, and collect lhs responses */
658: (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m);
659: PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge);
661: /* set lsh of hc, fire, and collect rhs responses */
662: (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m);
663: PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge);
665: /* count number of dofs I own that have signal and not in sep set */
666: for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) {
667: if (!used[i]) {
668: /* number of unmarked dofs on node */
669: ct++;
670: /* number of dofs to be marked on lhs hc */
671: if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++;
672: /* number of dofs to be marked on rhs hc */
673: if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++;
674: }
675: }
677: /* for the non-symmetric case we need separators of width 2 */
678: /* so take both sides */
679: (id<mask) ? (sum[2]=ct) : (sum[3]=ct);
680: PCTFS_giop_hc(sum,w,4,op,edge);
682: ct=0;
683: if (id<mask) {
684: /* mark dofs I own that have signal and not in sep set */
685: for (i=0;i<n;i++) {
686: if ((!used[i])&&(lhs[i]!=0.0)) {
687: ct++; nfo++;
688: *--iptr = local2global[i];
689: used[i] =edge;
690: }
691: }
692: /* LSH hc summation of ct should be sum[0] */
693: } else {
694: /* mark dofs I own that have signal and not in sep set */
695: for (i=0; i<n; i++) {
696: if ((!used[i])&&(rhs[i]!=0.0)) {
697: ct++; nfo++;
698: *--iptr = local2global[i];
699: used[i] = edge;
700: }
701: }
702: /* RSH hc summation of ct should be sum[1] */
703: }
705: if (ct>1) PCTFS_ivec_sort(iptr,ct);
706: lnsep[edge]=ct;
707: nsep[edge] =sum[0]+sum[1];
708: dir [edge] =BOTH;
710: /* LATER or we can recur on these to order seps at this level */
711: /* do we need full set of separators for this? */
713: /* fold rhs hc into lower */
714: if (id>=mask) id-=mask;
715: }
716: }
718: /* level 0 is on processor case - so mark the remainder */
719: for (ct=i=0;i<n;i++) {
720: if (!used[i]) {
721: ct++; nfo++;
722: *--iptr = local2global[i];
723: used[i] = edge;
724: }
725: }
726: if (ct>1) PCTFS_ivec_sort(iptr,ct);
727: lnsep[edge]=ct;
728: nsep [edge]=ct;
729: dir [edge]=BOTH;
731: xyt_handle->info->nsep = nsep;
732: xyt_handle->info->lnsep = lnsep;
733: xyt_handle->info->fo = fo;
734: xyt_handle->info->nfo = nfo;
736: free(dir);
737: free(lhs);
738: free(rhs);
739: free(used);
740: return 0;
741: }
743: /**************************************xyt.c***********************************/
744: static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data)
745: {
746: mv_info *mvi;
748: mvi = (mv_info*)malloc(sizeof(mv_info));
749: mvi->n = n;
750: mvi->m = m;
751: mvi->n_global = -1;
752: mvi->m_global = -1;
753: mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt));
755: PCTFS_ivec_copy(mvi->local2global,local2global,m);
756: mvi->local2global[m] = INT_MAX;
757: mvi->matvec = matvec;
758: mvi->grid_data = grid_data;
760: /* set xyt communication handle to perform restricted matvec */
761: mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes);
763: return(mvi);
764: }
766: /**************************************xyt.c***********************************/
767: static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u)
768: {
769: A->matvec((mv_info*)A->grid_data,v,u);
770: return 0;
771: }