1827bd09bSSatish Balay 2827bd09bSSatish Balay /*************************************xyt.c************************************ 3827bd09bSSatish Balay Module Name: xyt 4827bd09bSSatish Balay Module Info: 5827bd09bSSatish Balay 6827bd09bSSatish Balay author: Henry M. Tufo III 7827bd09bSSatish Balay e-mail: hmt@asci.uchicago.edu 8827bd09bSSatish Balay contact: 9827bd09bSSatish Balay +--------------------------------+--------------------------------+ 10827bd09bSSatish Balay |MCS Division - Building 221 |Department of Computer Science | 11827bd09bSSatish Balay |Argonne National Laboratory |Ryerson 152 | 12827bd09bSSatish Balay |9700 S. Cass Avenue |The University of Chicago | 13827bd09bSSatish Balay |Argonne, IL 60439 |Chicago, IL 60637 | 14827bd09bSSatish Balay |(630) 252-5354/5986 ph/fx |(773) 702-6019/8487 ph/fx | 15827bd09bSSatish Balay +--------------------------------+--------------------------------+ 16827bd09bSSatish Balay 17827bd09bSSatish Balay Last Modification: 3.20.01 18827bd09bSSatish Balay **************************************xyt.c***********************************/ 19c6db04a5SJed Brown #include <../src/ksp/pc/impls/tfs/tfs.h> 20827bd09bSSatish Balay 21827bd09bSSatish Balay #define LEFT -1 22827bd09bSSatish Balay #define RIGHT 1 23827bd09bSSatish Balay #define BOTH 0 24827bd09bSSatish Balay 25827bd09bSSatish Balay typedef struct xyt_solver_info { 2652f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 2752f87cdaSBarry Smith PetscInt nnz, max_nnz, msg_buf_sz; 2852f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo, *stages; 2952f87cdaSBarry Smith PetscInt *xcol_sz, *xcol_indices; 30a501084fSBarry Smith PetscScalar **xcol_vals, *x, *solve_uu, *solve_w; 3152f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 32a501084fSBarry Smith PetscScalar **ycol_vals, *y; 3352f87cdaSBarry Smith PetscInt nsolves; 34a501084fSBarry Smith PetscScalar tot_solve_time; 35827bd09bSSatish Balay } xyt_info; 36827bd09bSSatish Balay 37827bd09bSSatish Balay 38827bd09bSSatish Balay typedef struct matvec_info { 3952f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 4052f87cdaSBarry Smith PetscInt *local2global; 41ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 42a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*); 43827bd09bSSatish Balay void *grid_data; 44827bd09bSSatish Balay } mv_info; 45827bd09bSSatish Balay 46827bd09bSSatish Balay struct xyt_CDT { 4752f87cdaSBarry Smith PetscInt id; 4852f87cdaSBarry Smith PetscInt ns; 4952f87cdaSBarry Smith PetscInt level; 50827bd09bSSatish Balay xyt_info *info; 51827bd09bSSatish Balay mv_info *mvi; 52827bd09bSSatish Balay }; 53827bd09bSSatish Balay 5452f87cdaSBarry Smith static PetscInt n_xyt =0; 5552f87cdaSBarry Smith static PetscInt n_xyt_handles=0; 56827bd09bSSatish Balay 57827bd09bSSatish Balay /* prototypes */ 583fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs); 593fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle); 603fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle); 613fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 62*38dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle); 63*38dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle); 645c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data); 65827bd09bSSatish Balay 667b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 677b1ae94cSBarry Smith xyt_ADT XYT_new(void) 68827bd09bSSatish Balay { 69827bd09bSSatish Balay xyt_ADT xyt_handle; 70827bd09bSSatish Balay 71827bd09bSSatish Balay /* rolling count on n_xyt ... pot. problem here */ 72827bd09bSSatish Balay n_xyt_handles++; 73a501084fSBarry Smith xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT)); 74827bd09bSSatish Balay xyt_handle->id = ++n_xyt; 75827bd09bSSatish Balay xyt_handle->info = NULL; 76827bd09bSSatish Balay xyt_handle->mvi = NULL; 77827bd09bSSatish Balay 78827bd09bSSatish Balay return(xyt_handle); 79827bd09bSSatish Balay } 80827bd09bSSatish Balay 817b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 82*38dcbb09SBarry Smith PetscErrorCode XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */ 8352f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 8452f87cdaSBarry Smith PetscInt n, /* local num rows */ 8552f87cdaSBarry Smith PetscInt m, /* local num cols */ 865c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */ 871147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 88827bd09bSSatish Balay { 89827bd09bSSatish Balay 90b1c944f5SJed Brown PCTFS_comm_init(); 91827bd09bSSatish Balay check_handle(xyt_handle); 92827bd09bSSatish Balay 93827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 94b1c944f5SJed Brown 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); 95827bd09bSSatish Balay 96827bd09bSSatish Balay /* space for X info */ 97a501084fSBarry Smith xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info)); 98827bd09bSSatish Balay 99827bd09bSSatish Balay /* set up matvec handles */ 1005c8f6a95SKarl Rupp xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data); 101827bd09bSSatish Balay 102827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 103827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 104827bd09bSSatish Balay xyt_handle->ns=0; 105827bd09bSSatish Balay 106827bd09bSSatish Balay /* determine separators and generate firing order - NB xyt info set here */ 107827bd09bSSatish Balay det_separators(xyt_handle); 108827bd09bSSatish Balay 109827bd09bSSatish Balay return(do_xyt_factor(xyt_handle)); 110827bd09bSSatish Balay } 111827bd09bSSatish Balay 1127b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 113*38dcbb09SBarry Smith PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b) 114827bd09bSSatish Balay { 115b1c944f5SJed Brown PCTFS_comm_init(); 116827bd09bSSatish Balay check_handle(xyt_handle); 117827bd09bSSatish Balay 118827bd09bSSatish Balay /* need to copy b into x? */ 11959bc5b24SSatish Balay if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n); 120*38dcbb09SBarry Smith return do_xyt_solve(xyt_handle,x); 121827bd09bSSatish Balay } 122827bd09bSSatish Balay 1237b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 124*38dcbb09SBarry Smith PetscErrorCode XYT_free(xyt_ADT xyt_handle) 125827bd09bSSatish Balay { 126b1c944f5SJed Brown PCTFS_comm_init(); 127827bd09bSSatish Balay check_handle(xyt_handle); 128827bd09bSSatish Balay n_xyt_handles--; 129827bd09bSSatish Balay 130a501084fSBarry Smith free(xyt_handle->info->nsep); 131a501084fSBarry Smith free(xyt_handle->info->lnsep); 132a501084fSBarry Smith free(xyt_handle->info->fo); 133a501084fSBarry Smith free(xyt_handle->info->stages); 134a501084fSBarry Smith free(xyt_handle->info->solve_uu); 135a501084fSBarry Smith free(xyt_handle->info->solve_w); 136a501084fSBarry Smith free(xyt_handle->info->x); 137a501084fSBarry Smith free(xyt_handle->info->xcol_vals); 138a501084fSBarry Smith free(xyt_handle->info->xcol_sz); 139a501084fSBarry Smith free(xyt_handle->info->xcol_indices); 140a501084fSBarry Smith free(xyt_handle->info->y); 141a501084fSBarry Smith free(xyt_handle->info->ycol_vals); 142a501084fSBarry Smith free(xyt_handle->info->ycol_sz); 143a501084fSBarry Smith free(xyt_handle->info->ycol_indices); 144a501084fSBarry Smith free(xyt_handle->info); 145a501084fSBarry Smith free(xyt_handle->mvi->local2global); 146ca8e9878SJed Brown PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle); 147a501084fSBarry Smith free(xyt_handle->mvi); 148a501084fSBarry Smith free(xyt_handle); 149827bd09bSSatish Balay 150827bd09bSSatish Balay 151827bd09bSSatish Balay /* if the check fails we nuke */ 152a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 153827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 154827bd09bSSatish Balay return(0); 155827bd09bSSatish Balay } 156827bd09bSSatish Balay 1577b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 158*38dcbb09SBarry Smith PetscErrorCode XYT_stats(xyt_ADT xyt_handle) 159827bd09bSSatish Balay { 16052f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 16152f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 16252f87cdaSBarry Smith PetscInt vals[9], work[9]; 163a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 164827bd09bSSatish Balay 165b1c944f5SJed Brown PCTFS_comm_init(); 166827bd09bSSatish Balay check_handle(xyt_handle); 167827bd09bSSatish Balay 168827bd09bSSatish Balay /* if factorization not done there are no stats */ 1697d0a6c19SBarry Smith if (!xyt_handle->info||!xyt_handle->mvi) { 170b1c944f5SJed Brown if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n"); 171827bd09bSSatish Balay return 1; 172827bd09bSSatish Balay } 173827bd09bSSatish Balay 174827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz; 175827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n; 176827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz; 177b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 178827bd09bSSatish Balay 1792fa5cd67SKarl Rupp fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++; 180b1c944f5SJed Brown PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 181827bd09bSSatish Balay 182b1c944f5SJed Brown if (!PCTFS_my_id) { 183b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_nnz=%D\n",PCTFS_my_id,vals[0]); 184b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_nnz=%D\n",PCTFS_my_id,vals[1]); 185b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes); 186b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_nnz=%D\n",PCTFS_my_id,vals[2]); 18777b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(2d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5))); 18877b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(3d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667))); 189b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_n =%D\n",PCTFS_my_id,vals[3]); 190b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_n =%D\n",PCTFS_my_id,vals[4]); 191b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes); 192b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_n =%D\n",PCTFS_my_id,vals[5]); 193b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_buf=%D\n",PCTFS_my_id,vals[6]); 194b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_buf=%D\n",PCTFS_my_id,vals[7]); 195b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes); 196b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_slv=%g\n",PCTFS_my_id,fvals[0]); 197b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_slv=%g\n",PCTFS_my_id,fvals[1]); 198b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes); 199827bd09bSSatish Balay } 200827bd09bSSatish Balay 201827bd09bSSatish Balay return(0); 202827bd09bSSatish Balay } 203827bd09bSSatish Balay 204827bd09bSSatish Balay 205827bd09bSSatish Balay /*************************************xyt.c************************************ 206827bd09bSSatish Balay 207827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 208827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 209827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 210827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 211827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 212827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 213ca8e9878SJed Brown PCTFS_gs_init/gop). 214827bd09bSSatish Balay 215827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 216827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 217827bd09bSSatish Balay **************************************xyt.c***********************************/ 218*38dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle) 219827bd09bSSatish Balay { 2207b1ae94cSBarry Smith return xyt_generate(xyt_handle); 221827bd09bSSatish Balay } 222827bd09bSSatish Balay 2237b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 224*38dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle) 225827bd09bSSatish Balay { 22652f87cdaSBarry Smith PetscInt i,j,k,idx; 22752f87cdaSBarry Smith PetscInt dim, col; 228a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22952f87cdaSBarry Smith PetscInt *segs; 23052f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 23152f87cdaSBarry Smith PetscInt off, len; 232a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr; 23352f87cdaSBarry Smith PetscInt *iptr, flag; 23452f87cdaSBarry Smith PetscInt start =0, end, work; 23552f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 236ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 23752f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 23852f87cdaSBarry Smith PetscInt a_n =xyt_handle->mvi->n; 23952f87cdaSBarry Smith PetscInt a_m =xyt_handle->mvi->m; 24052f87cdaSBarry Smith PetscInt *a_local2global=xyt_handle->mvi->local2global; 24152f87cdaSBarry Smith PetscInt level; 24252f87cdaSBarry Smith PetscInt n, m; 24352f87cdaSBarry Smith PetscInt *xcol_sz, *xcol_indices, *stages; 244a501084fSBarry Smith PetscScalar **xcol_vals, *x; 24552f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 246a501084fSBarry Smith PetscScalar **ycol_vals, *y; 24752f87cdaSBarry Smith PetscInt n_global; 24852f87cdaSBarry Smith PetscInt xt_nnz =0, xt_max_nnz=0; 24952f87cdaSBarry Smith PetscInt yt_nnz =0, yt_max_nnz=0; 25052f87cdaSBarry Smith PetscInt xt_zero_nnz =0; 25152f87cdaSBarry Smith PetscInt xt_zero_nnz_0=0; 25252f87cdaSBarry Smith PetscInt yt_zero_nnz =0; 25352f87cdaSBarry Smith PetscInt yt_zero_nnz_0=0; 2544a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 255a501084fSBarry Smith PetscScalar dm1 = -1.0; 256d1528f56SBarry Smith PetscErrorCode ierr; 257827bd09bSSatish Balay 258827bd09bSSatish Balay n =xyt_handle->mvi->n; 259827bd09bSSatish Balay nsep =xyt_handle->info->nsep; 260827bd09bSSatish Balay lnsep =xyt_handle->info->lnsep; 261827bd09bSSatish Balay fo =xyt_handle->info->fo; 262827bd09bSSatish Balay end =lnsep[0]; 263827bd09bSSatish Balay level =xyt_handle->level; 264ca8e9878SJed Brown PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 265827bd09bSSatish Balay 266827bd09bSSatish Balay /* is there a null space? */ 267827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2682fa5cd67SKarl Rupp for (i=0, j=0; i<=level; i++) j+=nsep[i]; 269827bd09bSSatish Balay 270827bd09bSSatish Balay m = j-xyt_handle->ns; 27122d28d08SBarry Smith if (m!=j) { 27222d28d08SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);CHKERRQ(ierr); 27322d28d08SBarry Smith } 274827bd09bSSatish Balay 2754a0de3f6SBarry Smith ierr = PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);CHKERRQ(ierr); 276827bd09bSSatish Balay 277827bd09bSSatish Balay /* get and initialize storage for x local */ 278827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 27952f87cdaSBarry Smith xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 28052f87cdaSBarry Smith xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 281a501084fSBarry Smith xcol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 282db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 283827bd09bSSatish Balay xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1; 284827bd09bSSatish Balay xcol_vals[i] = NULL; 285827bd09bSSatish Balay } 286827bd09bSSatish Balay xcol_indices[j]=-1; 287827bd09bSSatish Balay 288827bd09bSSatish Balay /* get and initialize storage for y local */ 289827bd09bSSatish Balay /* note that y local is nxm and stored by columns */ 29052f87cdaSBarry Smith ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 29152f87cdaSBarry Smith ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 292a501084fSBarry Smith ycol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 293db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 294827bd09bSSatish Balay ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1; 295827bd09bSSatish Balay ycol_vals[i] = NULL; 296827bd09bSSatish Balay } 297827bd09bSSatish Balay ycol_indices[j]=-1; 298827bd09bSSatish Balay 299827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 300827bd09bSSatish Balay /* this looks like nsep[]=segments */ 30152f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 30252f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 303ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 304ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 3052fa5cd67SKarl Rupp for (i=0; i<level; i++) segs[i+1] += segs[i]; 306827bd09bSSatish Balay stages[0] = segs[0]; 307827bd09bSSatish Balay 308827bd09bSSatish Balay /* temporary vectors */ 309a501084fSBarry Smith u = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 310a501084fSBarry Smith z = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 311a501084fSBarry Smith v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar)); 312a501084fSBarry Smith uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 313a501084fSBarry Smith w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 314827bd09bSSatish Balay 315827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 3162fa5cd67SKarl Rupp for (i=1, j=0; i<=level; i++) j+=nsep[i]; 317827bd09bSSatish Balay 318827bd09bSSatish Balay /* storage for sparse x values */ 319827bd09bSSatish Balay n_global = xyt_handle->info->n_global; 32085ec1a3cSBarry Smith xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 321a501084fSBarry Smith x = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 322a501084fSBarry Smith y = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 323827bd09bSSatish Balay 324827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 325827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 326db4deed7SKarl Rupp for (dim=i=j=0; i<m; i++) { 327827bd09bSSatish Balay /* time to move to the next level? */ 328d4af0d30SBarry Smith while (i==segs[dim]) { 329e32f2f54SBarry Smith if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 330827bd09bSSatish Balay stages[dim++]=i; 331827bd09bSSatish Balay end +=lnsep[dim]; 332827bd09bSSatish Balay } 333827bd09bSSatish Balay stages[dim]=i; 334827bd09bSSatish Balay 335827bd09bSSatish Balay /* which column are we firing? */ 336827bd09bSSatish Balay /* i.e. set v_l */ 337827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 338827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 339b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 340827bd09bSSatish Balay 341827bd09bSSatish Balay /* shouldn't need this */ 342db4deed7SKarl Rupp if (col==INT_MAX) { 343f1ed62a8SBarry Smith ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 344827bd09bSSatish Balay continue; 345827bd09bSSatish Balay } 346827bd09bSSatish Balay 347827bd09bSSatish Balay /* do I own it? I should */ 348ca8e9878SJed Brown PCTFS_rvec_zero(v,a_m); 3492fa5cd67SKarl Rupp if (col==fo[start]) { 350827bd09bSSatish Balay start++; 351ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_n); 3522fa5cd67SKarl Rupp if (idx!=-1) { 3532fa5cd67SKarl Rupp v[idx] = 1.0; 3542fa5cd67SKarl Rupp j++; 3552fa5cd67SKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n"); 356db4deed7SKarl Rupp } else { 357ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_m); 3582fa5cd67SKarl Rupp if (idx!=-1) v[idx] = 1.0; 359827bd09bSSatish Balay } 360827bd09bSSatish Balay 361827bd09bSSatish Balay /* perform u = A.v_l */ 362ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 363827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 364827bd09bSSatish Balay 365827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 366827bd09bSSatish Balay /* technically only need to zero out first i entries */ 367827bd09bSSatish Balay /* later turn this into an XYT_solve call ? */ 368ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 369827bd09bSSatish Balay y_ptr=y; 370827bd09bSSatish Balay iptr = ycol_indices; 371db4deed7SKarl Rupp for (k=0; k<i; k++) { 372827bd09bSSatish Balay off = *iptr++; 373827bd09bSSatish Balay len = *iptr++; 374c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3758b83055fSJed Brown PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1)); 376827bd09bSSatish Balay y_ptr+=len; 377827bd09bSSatish Balay } 378827bd09bSSatish Balay 379827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 380a00740a5SBarry Smith ierr = PCTFS_ssgl_radd (uu, w, dim, stages);CHKERRQ(ierr); 381827bd09bSSatish Balay 382827bd09bSSatish Balay /* z = X.uu */ 383ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 384827bd09bSSatish Balay x_ptr=x; 385827bd09bSSatish Balay iptr = xcol_indices; 386db4deed7SKarl Rupp for (k=0; k<i; k++) { 387827bd09bSSatish Balay off = *iptr++; 388827bd09bSSatish Balay len = *iptr++; 389c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3908b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1)); 391827bd09bSSatish Balay x_ptr+=len; 392827bd09bSSatish Balay } 393827bd09bSSatish Balay 394827bd09bSSatish Balay /* compute v_l = v_l - z */ 395ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 396c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 3978b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1)); 398827bd09bSSatish Balay 399827bd09bSSatish Balay /* compute u_l = A.v_l */ 4002fa5cd67SKarl Rupp if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); 401ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 402827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 403827bd09bSSatish Balay 404827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 405c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 4068b83055fSJed Brown PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,u,&i1)); 407827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 408b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 409827bd09bSSatish Balay 4108f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 411827bd09bSSatish Balay 412827bd09bSSatish Balay /* check for small alpha */ 413827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 41477b4d14cSPeter Brune if (PetscAbsScalar(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 415827bd09bSSatish Balay 416827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 417ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 418ca8e9878SJed Brown PCTFS_rvec_scale(u,1.0/alpha,n); 419827bd09bSSatish Balay 420827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 421827bd09bSSatish Balay flag = 1; 422827bd09bSSatish Balay off =len=0; 423db4deed7SKarl Rupp for (k=0; k<n; k++) { 424db4deed7SKarl Rupp if (v[k]!=0.0) { 425827bd09bSSatish Balay len=k; 426db4deed7SKarl Rupp if (flag) {off=k; flag=0;} 427827bd09bSSatish Balay } 428827bd09bSSatish Balay } 429827bd09bSSatish Balay 430827bd09bSSatish Balay len -= (off-1); 431827bd09bSSatish Balay 432db4deed7SKarl Rupp if (len>0) { 433db4deed7SKarl Rupp if ((xt_nnz+len)>xt_max_nnz) { 434f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 435827bd09bSSatish Balay xt_max_nnz *= 2; 436a501084fSBarry Smith x_ptr = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 437ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xt_nnz); 438a501084fSBarry Smith free(x); 439827bd09bSSatish Balay x = x_ptr; 440827bd09bSSatish Balay x_ptr+=xt_nnz; 441827bd09bSSatish Balay } 442827bd09bSSatish Balay xt_nnz += len; 443ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 444827bd09bSSatish Balay 445827bd09bSSatish Balay /* keep track of number of zeros */ 446db4deed7SKarl Rupp if (dim) { 447db4deed7SKarl Rupp for (k=0; k<len; k++) { 4482fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz++; 449827bd09bSSatish Balay } 450db4deed7SKarl Rupp } else { 451db4deed7SKarl Rupp for (k=0; k<len; k++) { 4522fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz_0++; 453827bd09bSSatish Balay } 454827bd09bSSatish Balay } 455827bd09bSSatish Balay xcol_indices[2*i] = off; 456827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = len; 457827bd09bSSatish Balay xcol_vals[i] = x_ptr; 458db4deed7SKarl Rupp } else { 459827bd09bSSatish Balay xcol_indices[2*i] = 0; 460827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = 0; 461827bd09bSSatish Balay xcol_vals[i] = x_ptr; 462827bd09bSSatish Balay } 463827bd09bSSatish Balay 464827bd09bSSatish Balay 465827bd09bSSatish Balay /* add newly generated column, u_l, to Y */ 466827bd09bSSatish Balay flag = 1; 467827bd09bSSatish Balay off =len=0; 468db4deed7SKarl Rupp for (k=0; k<n; k++) { 469db4deed7SKarl Rupp if (u[k]!=0.0) { 470827bd09bSSatish Balay len=k; 471db4deed7SKarl Rupp if (flag) { off=k; flag=0; } 472827bd09bSSatish Balay } 473827bd09bSSatish Balay } 474827bd09bSSatish Balay 475827bd09bSSatish Balay len -= (off-1); 476827bd09bSSatish Balay 477db4deed7SKarl Rupp if (len>0) { 478db4deed7SKarl Rupp if ((yt_nnz+len)>yt_max_nnz) { 479f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr); 480827bd09bSSatish Balay yt_max_nnz *= 2; 481a501084fSBarry Smith y_ptr = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 482ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,y,yt_nnz); 483a501084fSBarry Smith free(y); 484827bd09bSSatish Balay y = y_ptr; 485827bd09bSSatish Balay y_ptr+=yt_nnz; 486827bd09bSSatish Balay } 487827bd09bSSatish Balay yt_nnz += len; 488ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,u+off,len); 489827bd09bSSatish Balay 490827bd09bSSatish Balay /* keep track of number of zeros */ 491db4deed7SKarl Rupp if (dim) { 492db4deed7SKarl Rupp for (k=0; k<len; k++) { 4932fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz++; 494827bd09bSSatish Balay } 495db4deed7SKarl Rupp } else { 496db4deed7SKarl Rupp for (k=0; k<len; k++) { 4972fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz_0++; 498827bd09bSSatish Balay } 499827bd09bSSatish Balay } 500827bd09bSSatish Balay ycol_indices[2*i] = off; 501827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = len; 502827bd09bSSatish Balay ycol_vals[i] = y_ptr; 503db4deed7SKarl Rupp } else { 504827bd09bSSatish Balay ycol_indices[2*i] = 0; 505827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = 0; 506827bd09bSSatish Balay ycol_vals[i] = y_ptr; 507827bd09bSSatish Balay } 508827bd09bSSatish Balay } 509827bd09bSSatish Balay 510827bd09bSSatish Balay /* close off stages for execution phase */ 511db4deed7SKarl Rupp while (dim!=level) { 512827bd09bSSatish Balay stages[dim++]=i; 5134a0de3f6SBarry Smith ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 514827bd09bSSatish Balay } 515827bd09bSSatish Balay stages[dim]=i; 516827bd09bSSatish Balay 517827bd09bSSatish Balay xyt_handle->info->n =xyt_handle->mvi->n; 518827bd09bSSatish Balay xyt_handle->info->m =m; 519827bd09bSSatish Balay xyt_handle->info->nnz =xt_nnz + yt_nnz; 520827bd09bSSatish Balay xyt_handle->info->max_nnz =xt_max_nnz + yt_max_nnz; 521827bd09bSSatish Balay xyt_handle->info->msg_buf_sz =stages[level]-stages[0]; 522a501084fSBarry Smith xyt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 523a501084fSBarry Smith xyt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 524827bd09bSSatish Balay xyt_handle->info->x =x; 525827bd09bSSatish Balay xyt_handle->info->xcol_vals =xcol_vals; 526827bd09bSSatish Balay xyt_handle->info->xcol_sz =xcol_sz; 527827bd09bSSatish Balay xyt_handle->info->xcol_indices=xcol_indices; 528827bd09bSSatish Balay xyt_handle->info->stages =stages; 529827bd09bSSatish Balay xyt_handle->info->y =y; 530827bd09bSSatish Balay xyt_handle->info->ycol_vals =ycol_vals; 531827bd09bSSatish Balay xyt_handle->info->ycol_sz =ycol_sz; 532827bd09bSSatish Balay xyt_handle->info->ycol_indices=ycol_indices; 533827bd09bSSatish Balay 534a501084fSBarry Smith free(segs); 535a501084fSBarry Smith free(u); 536a501084fSBarry Smith free(v); 537a501084fSBarry Smith free(uu); 538a501084fSBarry Smith free(z); 539a501084fSBarry Smith free(w); 540827bd09bSSatish Balay 541827bd09bSSatish Balay return(0); 542827bd09bSSatish Balay } 543827bd09bSSatish Balay 5447b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5450924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 546827bd09bSSatish Balay { 54752f87cdaSBarry Smith PetscInt off, len, *iptr; 54852f87cdaSBarry Smith PetscInt level =xyt_handle->level; 54952f87cdaSBarry Smith PetscInt n =xyt_handle->info->n; 55052f87cdaSBarry Smith PetscInt m =xyt_handle->info->m; 55152f87cdaSBarry Smith PetscInt *stages =xyt_handle->info->stages; 55252f87cdaSBarry Smith PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 55352f87cdaSBarry Smith PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 554a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr, *uu_ptr; 555a501084fSBarry Smith PetscScalar *solve_uu=xyt_handle->info->solve_uu; 556a501084fSBarry Smith PetscScalar *solve_w =xyt_handle->info->solve_w; 557a501084fSBarry Smith PetscScalar *x =xyt_handle->info->x; 558a501084fSBarry Smith PetscScalar *y =xyt_handle->info->y; 5594a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 560c5df96a5SBarry Smith PetscErrorCode ierr; 561827bd09bSSatish Balay 5620924e98cSBarry Smith PetscFunctionBegin; 563827bd09bSSatish Balay uu_ptr=solve_uu; 564ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 565827bd09bSSatish Balay 566827bd09bSSatish Balay /* x = X.Y^T.b */ 567827bd09bSSatish Balay /* uu = Y^T.b */ 568947b95d8SBarry Smith for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) { 569c5df96a5SBarry Smith off =*iptr++; 570c5df96a5SBarry Smith len =*iptr++; 571c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 5728b83055fSJed Brown PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1)); 573827bd09bSSatish Balay } 574827bd09bSSatish Balay 575827bd09bSSatish Balay /* comunication of beta */ 576827bd09bSSatish Balay uu_ptr=solve_uu; 577a00740a5SBarry Smith if (level) {ierr = PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);CHKERRQ(ierr);} 578ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 579827bd09bSSatish Balay 580827bd09bSSatish Balay /* x = X.uu */ 581db4deed7SKarl Rupp for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) { 582c5df96a5SBarry Smith off =*iptr++; 583c5df96a5SBarry Smith len =*iptr++; 584c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 5858b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1)); 586827bd09bSSatish Balay } 5870924e98cSBarry Smith PetscFunctionReturn(0); 588827bd09bSSatish Balay } 589827bd09bSSatish Balay 5907b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5910924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle) 592827bd09bSSatish Balay { 59352f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 594827bd09bSSatish Balay 5950924e98cSBarry Smith PetscFunctionBegin; 596e32f2f54SBarry Smith if (!xyt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle); 597827bd09bSSatish Balay 598827bd09bSSatish Balay vals[0]=vals[1]=xyt_handle->id; 599b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 600e32f2f54SBarry Smith 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); 6010924e98cSBarry Smith PetscFunctionReturn(0); 602827bd09bSSatish Balay } 603827bd09bSSatish Balay 6047b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 6050924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle) 606827bd09bSSatish Balay { 60752f87cdaSBarry Smith PetscInt i, ct, id; 60852f87cdaSBarry Smith PetscInt mask, edge, *iptr; 60952f87cdaSBarry Smith PetscInt *dir, *used; 61052f87cdaSBarry Smith PetscInt sum[4], w[4]; 611a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 61252f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 613a501084fSBarry Smith PetscScalar *lhs, *rhs; 61452f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 615ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 61652f87cdaSBarry Smith PetscInt *local2global =xyt_handle->mvi->local2global; 61752f87cdaSBarry Smith PetscInt n =xyt_handle->mvi->n; 61852f87cdaSBarry Smith PetscInt m =xyt_handle->mvi->m; 61952f87cdaSBarry Smith PetscInt level =xyt_handle->level; 620ab824b78SBarry Smith PetscInt shared =0; 621d1528f56SBarry Smith PetscErrorCode ierr; 622827bd09bSSatish Balay 6230924e98cSBarry Smith PetscFunctionBegin; 62452f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62552f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62652f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62752f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 62852f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 629827bd09bSSatish Balay 630ca8e9878SJed Brown PCTFS_ivec_zero(dir,level+1); 631ca8e9878SJed Brown PCTFS_ivec_zero(nsep,level+1); 632ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 633ca8e9878SJed Brown PCTFS_ivec_set (fo,-1,n+1); 634ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 635827bd09bSSatish Balay 6368cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 6378cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 638827bd09bSSatish Balay 639827bd09bSSatish Balay /* determine the # of unique dof */ 640ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 641ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 642ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 643ca8e9878SJed Brown ierr = PetscInfo(0,"done first PCTFS_gs_gop_hc\n");CHKERRQ(ierr); 644ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 645947b95d8SBarry Smith for (i=0; i<n; i++) { 646db4deed7SKarl Rupp if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; } 6472fa5cd67SKarl Rupp if (lhs[i]!=1.0) shared=1; 648827bd09bSSatish Balay } 649827bd09bSSatish Balay 650b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 651827bd09bSSatish Balay rsum[0]+=0.1; 652827bd09bSSatish Balay rsum[1]+=0.1; 653827bd09bSSatish Balay 65452f87cdaSBarry Smith xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 65552f87cdaSBarry Smith xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 656827bd09bSSatish Balay 657827bd09bSSatish Balay /* determine separator sets top down */ 6582fa5cd67SKarl Rupp if (shared) { 659827bd09bSSatish Balay /* solution is to do as in the symmetric shared case but then */ 660827bd09bSSatish Balay /* pick the sub-hc with the most free dofs and do a mat-vec */ 661827bd09bSSatish Balay /* and pick up the responses on the other sub-hc from the */ 662827bd09bSSatish Balay /* initial separator set obtained from the symm. shared case */ 663e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n"); 6644e619c20SJed Brown /* [dead code deleted since it is unlikely to be completed] */ 665db4deed7SKarl Rupp } else { 666db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 667827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 668ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 669ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 670827bd09bSSatish Balay 671827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 672ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 673ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 674827bd09bSSatish Balay 675827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 676db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 677db4deed7SKarl Rupp if (!used[i]) { 678827bd09bSSatish Balay /* number of unmarked dofs on node */ 679827bd09bSSatish Balay ct++; 680827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6812fa5cd67SKarl Rupp if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++; 682827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6832fa5cd67SKarl Rupp if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++; 684827bd09bSSatish Balay } 685827bd09bSSatish Balay } 686827bd09bSSatish Balay 687827bd09bSSatish Balay /* for the non-symmetric case we need separators of width 2 */ 688827bd09bSSatish Balay /* so take both sides */ 689827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 690b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 691827bd09bSSatish Balay 692827bd09bSSatish Balay ct=0; 693db4deed7SKarl Rupp if (id<mask) { 694827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 695db4deed7SKarl Rupp for (i=0;i<n;i++) { 696db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 697827bd09bSSatish Balay ct++; nfo++; 698827bd09bSSatish Balay *--iptr = local2global[i]; 699827bd09bSSatish Balay used[i] =edge; 700827bd09bSSatish Balay } 701827bd09bSSatish Balay } 702827bd09bSSatish Balay /* LSH hc summation of ct should be sum[0] */ 703db4deed7SKarl Rupp } else { 704827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 705db4deed7SKarl Rupp for (i=0; i<n; i++) { 706db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 707827bd09bSSatish Balay ct++; nfo++; 708827bd09bSSatish Balay *--iptr = local2global[i]; 709827bd09bSSatish Balay used[i] = edge; 710827bd09bSSatish Balay } 711827bd09bSSatish Balay } 712827bd09bSSatish Balay /* RSH hc summation of ct should be sum[1] */ 713827bd09bSSatish Balay } 714827bd09bSSatish Balay 7152fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 716827bd09bSSatish Balay lnsep[edge]=ct; 717827bd09bSSatish Balay nsep[edge] =sum[0]+sum[1]; 718827bd09bSSatish Balay dir [edge] =BOTH; 719827bd09bSSatish Balay 720827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 721827bd09bSSatish Balay /* do we need full set of separators for this? */ 722827bd09bSSatish Balay 723827bd09bSSatish Balay /* fold rhs hc into lower */ 7242fa5cd67SKarl Rupp if (id>=mask) id-=mask; 725827bd09bSSatish Balay } 726827bd09bSSatish Balay } 727827bd09bSSatish Balay 728827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 729db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 730db4deed7SKarl Rupp if (!used[i]) { 731827bd09bSSatish Balay ct++; nfo++; 732827bd09bSSatish Balay *--iptr = local2global[i]; 733827bd09bSSatish Balay used[i] = edge; 734827bd09bSSatish Balay } 735827bd09bSSatish Balay } 7362fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 737827bd09bSSatish Balay lnsep[edge]=ct; 738827bd09bSSatish Balay nsep [edge]=ct; 739827bd09bSSatish Balay dir [edge]=BOTH; 740827bd09bSSatish Balay 741827bd09bSSatish Balay xyt_handle->info->nsep = nsep; 742827bd09bSSatish Balay xyt_handle->info->lnsep = lnsep; 743827bd09bSSatish Balay xyt_handle->info->fo = fo; 744827bd09bSSatish Balay xyt_handle->info->nfo = nfo; 745827bd09bSSatish Balay 746a501084fSBarry Smith free(dir); 747a501084fSBarry Smith free(lhs); 748a501084fSBarry Smith free(rhs); 749a501084fSBarry Smith free(used); 7500924e98cSBarry Smith PetscFunctionReturn(0); 751827bd09bSSatish Balay } 752827bd09bSSatish Balay 7537b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7545c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data) 755827bd09bSSatish Balay { 756827bd09bSSatish Balay mv_info *mvi; 757827bd09bSSatish Balay 758a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 759827bd09bSSatish Balay mvi->n = n; 760827bd09bSSatish Balay mvi->m = m; 761827bd09bSSatish Balay mvi->n_global = -1; 762827bd09bSSatish Balay mvi->m_global = -1; 76352f87cdaSBarry Smith mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt)); 7642fa5cd67SKarl Rupp 765ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 766827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7675c8f6a95SKarl Rupp mvi->matvec = matvec; 768827bd09bSSatish Balay mvi->grid_data = grid_data; 769827bd09bSSatish Balay 770827bd09bSSatish Balay /* set xyt communication handle to perform restricted matvec */ 771ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 772827bd09bSSatish Balay 773827bd09bSSatish Balay return(mvi); 774827bd09bSSatish Balay } 775827bd09bSSatish Balay 7767b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7773fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 778827bd09bSSatish Balay { 7790924e98cSBarry Smith PetscFunctionBegin; 780827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 7810924e98cSBarry Smith PetscFunctionReturn(0); 782827bd09bSSatish Balay } 783827bd09bSSatish Balay 784827bd09bSSatish Balay 785827bd09bSSatish Balay 786