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