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); 6252f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle); 6352f87cdaSBarry Smith static PetscInt 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***********************************/ 8252f87cdaSBarry Smith PetscInt 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***********************************/ 1138cda6cd7SBarry Smith PetscInt 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? */ 119*59bc5b24SSatish Balay if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n); 120*59bc5b24SSatish Balay do_xyt_solve(xyt_handle,x); 121827bd09bSSatish Balay 122827bd09bSSatish Balay return(0); 123827bd09bSSatish Balay } 124827bd09bSSatish Balay 1257b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 12652f87cdaSBarry Smith PetscInt XYT_free(xyt_ADT xyt_handle) 127827bd09bSSatish Balay { 128b1c944f5SJed Brown PCTFS_comm_init(); 129827bd09bSSatish Balay check_handle(xyt_handle); 130827bd09bSSatish Balay n_xyt_handles--; 131827bd09bSSatish Balay 132a501084fSBarry Smith free(xyt_handle->info->nsep); 133a501084fSBarry Smith free(xyt_handle->info->lnsep); 134a501084fSBarry Smith free(xyt_handle->info->fo); 135a501084fSBarry Smith free(xyt_handle->info->stages); 136a501084fSBarry Smith free(xyt_handle->info->solve_uu); 137a501084fSBarry Smith free(xyt_handle->info->solve_w); 138a501084fSBarry Smith free(xyt_handle->info->x); 139a501084fSBarry Smith free(xyt_handle->info->xcol_vals); 140a501084fSBarry Smith free(xyt_handle->info->xcol_sz); 141a501084fSBarry Smith free(xyt_handle->info->xcol_indices); 142a501084fSBarry Smith free(xyt_handle->info->y); 143a501084fSBarry Smith free(xyt_handle->info->ycol_vals); 144a501084fSBarry Smith free(xyt_handle->info->ycol_sz); 145a501084fSBarry Smith free(xyt_handle->info->ycol_indices); 146a501084fSBarry Smith free(xyt_handle->info); 147a501084fSBarry Smith free(xyt_handle->mvi->local2global); 148ca8e9878SJed Brown PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle); 149a501084fSBarry Smith free(xyt_handle->mvi); 150a501084fSBarry Smith free(xyt_handle); 151827bd09bSSatish Balay 152827bd09bSSatish Balay 153827bd09bSSatish Balay /* if the check fails we nuke */ 154a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 155827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 156827bd09bSSatish Balay return(0); 157827bd09bSSatish Balay } 158827bd09bSSatish Balay 1597b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 16052f87cdaSBarry Smith PetscInt XYT_stats(xyt_ADT xyt_handle) 161827bd09bSSatish Balay { 16252f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 16352f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 16452f87cdaSBarry Smith PetscInt vals[9], work[9]; 165a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 166827bd09bSSatish Balay 167b1c944f5SJed Brown PCTFS_comm_init(); 168827bd09bSSatish Balay check_handle(xyt_handle); 169827bd09bSSatish Balay 170827bd09bSSatish Balay /* if factorization not done there are no stats */ 1717d0a6c19SBarry Smith if (!xyt_handle->info||!xyt_handle->mvi) { 172b1c944f5SJed Brown if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n"); 173827bd09bSSatish Balay return 1; 174827bd09bSSatish Balay } 175827bd09bSSatish Balay 176827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz; 177827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n; 178827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz; 179b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 180827bd09bSSatish Balay 1812fa5cd67SKarl Rupp fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++; 182b1c944f5SJed Brown PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 183827bd09bSSatish Balay 184b1c944f5SJed Brown if (!PCTFS_my_id) { 185b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_nnz=%D\n",PCTFS_my_id,vals[0]); 186b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_nnz=%D\n",PCTFS_my_id,vals[1]); 187b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes); 188b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_nnz=%D\n",PCTFS_my_id,vals[2]); 18977b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(2d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5))); 19077b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(3d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667))); 191b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_n =%D\n",PCTFS_my_id,vals[3]); 192b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_n =%D\n",PCTFS_my_id,vals[4]); 193b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes); 194b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_n =%D\n",PCTFS_my_id,vals[5]); 195b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_buf=%D\n",PCTFS_my_id,vals[6]); 196b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_buf=%D\n",PCTFS_my_id,vals[7]); 197b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes); 198b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_slv=%g\n",PCTFS_my_id,fvals[0]); 199b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_slv=%g\n",PCTFS_my_id,fvals[1]); 200b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes); 201827bd09bSSatish Balay } 202827bd09bSSatish Balay 203827bd09bSSatish Balay return(0); 204827bd09bSSatish Balay } 205827bd09bSSatish Balay 206827bd09bSSatish Balay 207827bd09bSSatish Balay /*************************************xyt.c************************************ 208827bd09bSSatish Balay 209827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 210827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 211827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 212827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 213827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 214827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 215ca8e9878SJed Brown PCTFS_gs_init/gop). 216827bd09bSSatish Balay 217827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 218827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 219827bd09bSSatish Balay **************************************xyt.c***********************************/ 22052f87cdaSBarry Smith static PetscInt do_xyt_factor(xyt_ADT xyt_handle) 221827bd09bSSatish Balay { 2227b1ae94cSBarry Smith return xyt_generate(xyt_handle); 223827bd09bSSatish Balay } 224827bd09bSSatish Balay 2257b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 22652f87cdaSBarry Smith static PetscInt xyt_generate(xyt_ADT xyt_handle) 227827bd09bSSatish Balay { 22852f87cdaSBarry Smith PetscInt i,j,k,idx; 22952f87cdaSBarry Smith PetscInt dim, col; 230a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 23152f87cdaSBarry Smith PetscInt *segs; 23252f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 23352f87cdaSBarry Smith PetscInt off, len; 234a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr; 23552f87cdaSBarry Smith PetscInt *iptr, flag; 23652f87cdaSBarry Smith PetscInt start =0, end, work; 23752f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 238ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 23952f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 24052f87cdaSBarry Smith PetscInt a_n =xyt_handle->mvi->n; 24152f87cdaSBarry Smith PetscInt a_m =xyt_handle->mvi->m; 24252f87cdaSBarry Smith PetscInt *a_local2global=xyt_handle->mvi->local2global; 24352f87cdaSBarry Smith PetscInt level; 24452f87cdaSBarry Smith PetscInt n, m; 24552f87cdaSBarry Smith PetscInt *xcol_sz, *xcol_indices, *stages; 246a501084fSBarry Smith PetscScalar **xcol_vals, *x; 24752f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 248a501084fSBarry Smith PetscScalar **ycol_vals, *y; 24952f87cdaSBarry Smith PetscInt n_global; 25052f87cdaSBarry Smith PetscInt xt_nnz =0, xt_max_nnz=0; 25152f87cdaSBarry Smith PetscInt yt_nnz =0, yt_max_nnz=0; 25252f87cdaSBarry Smith PetscInt xt_zero_nnz =0; 25352f87cdaSBarry Smith PetscInt xt_zero_nnz_0=0; 25452f87cdaSBarry Smith PetscInt yt_zero_nnz =0; 25552f87cdaSBarry Smith PetscInt yt_zero_nnz_0=0; 2564a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 257a501084fSBarry Smith PetscScalar dm1 = -1.0; 258d1528f56SBarry Smith PetscErrorCode ierr; 259827bd09bSSatish Balay 260827bd09bSSatish Balay n =xyt_handle->mvi->n; 261827bd09bSSatish Balay nsep =xyt_handle->info->nsep; 262827bd09bSSatish Balay lnsep =xyt_handle->info->lnsep; 263827bd09bSSatish Balay fo =xyt_handle->info->fo; 264827bd09bSSatish Balay end =lnsep[0]; 265827bd09bSSatish Balay level =xyt_handle->level; 266ca8e9878SJed Brown PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 267827bd09bSSatish Balay 268827bd09bSSatish Balay /* is there a null space? */ 269827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2702fa5cd67SKarl Rupp for (i=0, j=0; i<=level; i++) j+=nsep[i]; 271827bd09bSSatish Balay 272827bd09bSSatish Balay m = j-xyt_handle->ns; 27322d28d08SBarry Smith if (m!=j) { 27422d28d08SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns);CHKERRQ(ierr); 27522d28d08SBarry Smith } 276827bd09bSSatish Balay 2774a0de3f6SBarry Smith ierr = PetscInfo2(0,"xyt_generate() :: X(%D,%D)\n",n,m);CHKERRQ(ierr); 278827bd09bSSatish Balay 279827bd09bSSatish Balay /* get and initialize storage for x local */ 280827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 28152f87cdaSBarry Smith xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 28252f87cdaSBarry Smith xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 283a501084fSBarry Smith xcol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 284db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 285827bd09bSSatish Balay xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1; 286827bd09bSSatish Balay xcol_vals[i] = NULL; 287827bd09bSSatish Balay } 288827bd09bSSatish Balay xcol_indices[j]=-1; 289827bd09bSSatish Balay 290827bd09bSSatish Balay /* get and initialize storage for y local */ 291827bd09bSSatish Balay /* note that y local is nxm and stored by columns */ 29252f87cdaSBarry Smith ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 29352f87cdaSBarry Smith ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 294a501084fSBarry Smith ycol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 295db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 296827bd09bSSatish Balay ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1; 297827bd09bSSatish Balay ycol_vals[i] = NULL; 298827bd09bSSatish Balay } 299827bd09bSSatish Balay ycol_indices[j]=-1; 300827bd09bSSatish Balay 301827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 302827bd09bSSatish Balay /* this looks like nsep[]=segments */ 30352f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 30452f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 305ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 306ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 3072fa5cd67SKarl Rupp for (i=0; i<level; i++) segs[i+1] += segs[i]; 308827bd09bSSatish Balay stages[0] = segs[0]; 309827bd09bSSatish Balay 310827bd09bSSatish Balay /* temporary vectors */ 311a501084fSBarry Smith u = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 312a501084fSBarry Smith z = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 313a501084fSBarry Smith v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar)); 314a501084fSBarry Smith uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 315a501084fSBarry Smith w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 316827bd09bSSatish Balay 317827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 3182fa5cd67SKarl Rupp for (i=1, j=0; i<=level; i++) j+=nsep[i]; 319827bd09bSSatish Balay 320827bd09bSSatish Balay /* storage for sparse x values */ 321827bd09bSSatish Balay n_global = xyt_handle->info->n_global; 32285ec1a3cSBarry Smith xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 323a501084fSBarry Smith x = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 324a501084fSBarry Smith y = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 325827bd09bSSatish Balay 326827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 327827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 328db4deed7SKarl Rupp for (dim=i=j=0; i<m; i++) { 329827bd09bSSatish Balay /* time to move to the next level? */ 330d4af0d30SBarry Smith while (i==segs[dim]) { 331e32f2f54SBarry Smith if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 332827bd09bSSatish Balay stages[dim++]=i; 333827bd09bSSatish Balay end +=lnsep[dim]; 334827bd09bSSatish Balay } 335827bd09bSSatish Balay stages[dim]=i; 336827bd09bSSatish Balay 337827bd09bSSatish Balay /* which column are we firing? */ 338827bd09bSSatish Balay /* i.e. set v_l */ 339827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 340827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 341b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 342827bd09bSSatish Balay 343827bd09bSSatish Balay /* shouldn't need this */ 344db4deed7SKarl Rupp if (col==INT_MAX) { 345f1ed62a8SBarry Smith ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 346827bd09bSSatish Balay continue; 347827bd09bSSatish Balay } 348827bd09bSSatish Balay 349827bd09bSSatish Balay /* do I own it? I should */ 350ca8e9878SJed Brown PCTFS_rvec_zero(v,a_m); 3512fa5cd67SKarl Rupp if (col==fo[start]) { 352827bd09bSSatish Balay start++; 353ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_n); 3542fa5cd67SKarl Rupp if (idx!=-1) { 3552fa5cd67SKarl Rupp v[idx] = 1.0; 3562fa5cd67SKarl Rupp j++; 3572fa5cd67SKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n"); 358db4deed7SKarl Rupp } else { 359ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_m); 3602fa5cd67SKarl Rupp if (idx!=-1) v[idx] = 1.0; 361827bd09bSSatish Balay } 362827bd09bSSatish Balay 363827bd09bSSatish Balay /* perform u = A.v_l */ 364ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 365827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 366827bd09bSSatish Balay 367827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 368827bd09bSSatish Balay /* technically only need to zero out first i entries */ 369827bd09bSSatish Balay /* later turn this into an XYT_solve call ? */ 370ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 371827bd09bSSatish Balay y_ptr=y; 372827bd09bSSatish Balay iptr = ycol_indices; 373db4deed7SKarl Rupp for (k=0; k<i; k++) { 374827bd09bSSatish Balay off = *iptr++; 375827bd09bSSatish Balay len = *iptr++; 376c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3778b83055fSJed Brown PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1)); 378827bd09bSSatish Balay y_ptr+=len; 379827bd09bSSatish Balay } 380827bd09bSSatish Balay 381827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 382b1c944f5SJed Brown PCTFS_ssgl_radd (uu, w, dim, stages); 383827bd09bSSatish Balay 384827bd09bSSatish Balay /* z = X.uu */ 385ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 386827bd09bSSatish Balay x_ptr=x; 387827bd09bSSatish Balay iptr = xcol_indices; 388db4deed7SKarl Rupp for (k=0; k<i; k++) { 389827bd09bSSatish Balay off = *iptr++; 390827bd09bSSatish Balay len = *iptr++; 391c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3928b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1)); 393827bd09bSSatish Balay x_ptr+=len; 394827bd09bSSatish Balay } 395827bd09bSSatish Balay 396827bd09bSSatish Balay /* compute v_l = v_l - z */ 397ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 398c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 3998b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1)); 400827bd09bSSatish Balay 401827bd09bSSatish Balay /* compute u_l = A.v_l */ 4022fa5cd67SKarl Rupp if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); 403ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 404827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 405827bd09bSSatish Balay 406827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 407c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 4088b83055fSJed Brown PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,u,&i1)); 409827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 410b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 411827bd09bSSatish Balay 4128f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 413827bd09bSSatish Balay 414827bd09bSSatish Balay /* check for small alpha */ 415827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 41677b4d14cSPeter Brune if (PetscAbsScalar(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 417827bd09bSSatish Balay 418827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 419ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 420ca8e9878SJed Brown PCTFS_rvec_scale(u,1.0/alpha,n); 421827bd09bSSatish Balay 422827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 423827bd09bSSatish Balay flag = 1; 424827bd09bSSatish Balay off =len=0; 425db4deed7SKarl Rupp for (k=0; k<n; k++) { 426db4deed7SKarl Rupp if (v[k]!=0.0) { 427827bd09bSSatish Balay len=k; 428db4deed7SKarl Rupp if (flag) {off=k; flag=0;} 429827bd09bSSatish Balay } 430827bd09bSSatish Balay } 431827bd09bSSatish Balay 432827bd09bSSatish Balay len -= (off-1); 433827bd09bSSatish Balay 434db4deed7SKarl Rupp if (len>0) { 435db4deed7SKarl Rupp if ((xt_nnz+len)>xt_max_nnz) { 436f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 437827bd09bSSatish Balay xt_max_nnz *= 2; 438a501084fSBarry Smith x_ptr = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 439ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xt_nnz); 440a501084fSBarry Smith free(x); 441827bd09bSSatish Balay x = x_ptr; 442827bd09bSSatish Balay x_ptr+=xt_nnz; 443827bd09bSSatish Balay } 444827bd09bSSatish Balay xt_nnz += len; 445ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 446827bd09bSSatish Balay 447827bd09bSSatish Balay /* keep track of number of zeros */ 448db4deed7SKarl Rupp if (dim) { 449db4deed7SKarl Rupp for (k=0; k<len; k++) { 4502fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz++; 451827bd09bSSatish Balay } 452db4deed7SKarl Rupp } else { 453db4deed7SKarl Rupp for (k=0; k<len; k++) { 4542fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz_0++; 455827bd09bSSatish Balay } 456827bd09bSSatish Balay } 457827bd09bSSatish Balay xcol_indices[2*i] = off; 458827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = len; 459827bd09bSSatish Balay xcol_vals[i] = x_ptr; 460db4deed7SKarl Rupp } else { 461827bd09bSSatish Balay xcol_indices[2*i] = 0; 462827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = 0; 463827bd09bSSatish Balay xcol_vals[i] = x_ptr; 464827bd09bSSatish Balay } 465827bd09bSSatish Balay 466827bd09bSSatish Balay 467827bd09bSSatish Balay /* add newly generated column, u_l, to Y */ 468827bd09bSSatish Balay flag = 1; 469827bd09bSSatish Balay off =len=0; 470db4deed7SKarl Rupp for (k=0; k<n; k++) { 471db4deed7SKarl Rupp if (u[k]!=0.0) { 472827bd09bSSatish Balay len=k; 473db4deed7SKarl Rupp if (flag) { off=k; flag=0; } 474827bd09bSSatish Balay } 475827bd09bSSatish Balay } 476827bd09bSSatish Balay 477827bd09bSSatish Balay len -= (off-1); 478827bd09bSSatish Balay 479db4deed7SKarl Rupp if (len>0) { 480db4deed7SKarl Rupp if ((yt_nnz+len)>yt_max_nnz) { 481f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for Y by 2x!\n");CHKERRQ(ierr); 482827bd09bSSatish Balay yt_max_nnz *= 2; 483a501084fSBarry Smith y_ptr = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 484ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,y,yt_nnz); 485a501084fSBarry Smith free(y); 486827bd09bSSatish Balay y = y_ptr; 487827bd09bSSatish Balay y_ptr+=yt_nnz; 488827bd09bSSatish Balay } 489827bd09bSSatish Balay yt_nnz += len; 490ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,u+off,len); 491827bd09bSSatish Balay 492827bd09bSSatish Balay /* keep track of number of zeros */ 493db4deed7SKarl Rupp if (dim) { 494db4deed7SKarl Rupp for (k=0; k<len; k++) { 4952fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz++; 496827bd09bSSatish Balay } 497db4deed7SKarl Rupp } else { 498db4deed7SKarl Rupp for (k=0; k<len; k++) { 4992fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz_0++; 500827bd09bSSatish Balay } 501827bd09bSSatish Balay } 502827bd09bSSatish Balay ycol_indices[2*i] = off; 503827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = len; 504827bd09bSSatish Balay ycol_vals[i] = y_ptr; 505db4deed7SKarl Rupp } else { 506827bd09bSSatish Balay ycol_indices[2*i] = 0; 507827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = 0; 508827bd09bSSatish Balay ycol_vals[i] = y_ptr; 509827bd09bSSatish Balay } 510827bd09bSSatish Balay } 511827bd09bSSatish Balay 512827bd09bSSatish Balay /* close off stages for execution phase */ 513db4deed7SKarl Rupp while (dim!=level) { 514827bd09bSSatish Balay stages[dim++]=i; 5154a0de3f6SBarry Smith ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 516827bd09bSSatish Balay } 517827bd09bSSatish Balay stages[dim]=i; 518827bd09bSSatish Balay 519827bd09bSSatish Balay xyt_handle->info->n =xyt_handle->mvi->n; 520827bd09bSSatish Balay xyt_handle->info->m =m; 521827bd09bSSatish Balay xyt_handle->info->nnz =xt_nnz + yt_nnz; 522827bd09bSSatish Balay xyt_handle->info->max_nnz =xt_max_nnz + yt_max_nnz; 523827bd09bSSatish Balay xyt_handle->info->msg_buf_sz =stages[level]-stages[0]; 524a501084fSBarry Smith xyt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 525a501084fSBarry Smith xyt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 526827bd09bSSatish Balay xyt_handle->info->x =x; 527827bd09bSSatish Balay xyt_handle->info->xcol_vals =xcol_vals; 528827bd09bSSatish Balay xyt_handle->info->xcol_sz =xcol_sz; 529827bd09bSSatish Balay xyt_handle->info->xcol_indices=xcol_indices; 530827bd09bSSatish Balay xyt_handle->info->stages =stages; 531827bd09bSSatish Balay xyt_handle->info->y =y; 532827bd09bSSatish Balay xyt_handle->info->ycol_vals =ycol_vals; 533827bd09bSSatish Balay xyt_handle->info->ycol_sz =ycol_sz; 534827bd09bSSatish Balay xyt_handle->info->ycol_indices=ycol_indices; 535827bd09bSSatish Balay 536a501084fSBarry Smith free(segs); 537a501084fSBarry Smith free(u); 538a501084fSBarry Smith free(v); 539a501084fSBarry Smith free(uu); 540a501084fSBarry Smith free(z); 541a501084fSBarry Smith free(w); 542827bd09bSSatish Balay 543827bd09bSSatish Balay return(0); 544827bd09bSSatish Balay } 545827bd09bSSatish Balay 5467b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5470924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 548827bd09bSSatish Balay { 54952f87cdaSBarry Smith PetscInt off, len, *iptr; 55052f87cdaSBarry Smith PetscInt level =xyt_handle->level; 55152f87cdaSBarry Smith PetscInt n =xyt_handle->info->n; 55252f87cdaSBarry Smith PetscInt m =xyt_handle->info->m; 55352f87cdaSBarry Smith PetscInt *stages =xyt_handle->info->stages; 55452f87cdaSBarry Smith PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 55552f87cdaSBarry Smith PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 556a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr, *uu_ptr; 557a501084fSBarry Smith PetscScalar *solve_uu=xyt_handle->info->solve_uu; 558a501084fSBarry Smith PetscScalar *solve_w =xyt_handle->info->solve_w; 559a501084fSBarry Smith PetscScalar *x =xyt_handle->info->x; 560a501084fSBarry Smith PetscScalar *y =xyt_handle->info->y; 5614a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 562c5df96a5SBarry Smith PetscErrorCode ierr; 563827bd09bSSatish Balay 5640924e98cSBarry Smith PetscFunctionBegin; 565827bd09bSSatish Balay uu_ptr=solve_uu; 566ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 567827bd09bSSatish Balay 568827bd09bSSatish Balay /* x = X.Y^T.b */ 569827bd09bSSatish Balay /* uu = Y^T.b */ 570947b95d8SBarry Smith for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) { 571c5df96a5SBarry Smith off =*iptr++; 572c5df96a5SBarry Smith len =*iptr++; 573c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 5748b83055fSJed Brown PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1)); 575827bd09bSSatish Balay } 576827bd09bSSatish Balay 577827bd09bSSatish Balay /* comunication of beta */ 578827bd09bSSatish Balay uu_ptr=solve_uu; 5792fa5cd67SKarl Rupp if (level) PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages); 580ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 581827bd09bSSatish Balay 582827bd09bSSatish Balay /* x = X.uu */ 583db4deed7SKarl Rupp for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) { 584c5df96a5SBarry Smith off =*iptr++; 585c5df96a5SBarry Smith len =*iptr++; 586c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 5878b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1)); 588827bd09bSSatish Balay } 5890924e98cSBarry Smith PetscFunctionReturn(0); 590827bd09bSSatish Balay } 591827bd09bSSatish Balay 5927b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5930924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle) 594827bd09bSSatish Balay { 59552f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 596827bd09bSSatish Balay 5970924e98cSBarry Smith PetscFunctionBegin; 598e32f2f54SBarry Smith if (!xyt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xyt_handle); 599827bd09bSSatish Balay 600827bd09bSSatish Balay vals[0]=vals[1]=xyt_handle->id; 601b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 602e32f2f54SBarry 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); 6030924e98cSBarry Smith PetscFunctionReturn(0); 604827bd09bSSatish Balay } 605827bd09bSSatish Balay 6067b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 6070924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle) 608827bd09bSSatish Balay { 60952f87cdaSBarry Smith PetscInt i, ct, id; 61052f87cdaSBarry Smith PetscInt mask, edge, *iptr; 61152f87cdaSBarry Smith PetscInt *dir, *used; 61252f87cdaSBarry Smith PetscInt sum[4], w[4]; 613a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 61452f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 615a501084fSBarry Smith PetscScalar *lhs, *rhs; 61652f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 617ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 61852f87cdaSBarry Smith PetscInt *local2global =xyt_handle->mvi->local2global; 61952f87cdaSBarry Smith PetscInt n =xyt_handle->mvi->n; 62052f87cdaSBarry Smith PetscInt m =xyt_handle->mvi->m; 62152f87cdaSBarry Smith PetscInt level =xyt_handle->level; 622ab824b78SBarry Smith PetscInt shared =0; 623d1528f56SBarry Smith PetscErrorCode ierr; 624827bd09bSSatish Balay 6250924e98cSBarry Smith PetscFunctionBegin; 62652f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62752f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62852f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62952f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 63052f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 631827bd09bSSatish Balay 632ca8e9878SJed Brown PCTFS_ivec_zero(dir,level+1); 633ca8e9878SJed Brown PCTFS_ivec_zero(nsep,level+1); 634ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 635ca8e9878SJed Brown PCTFS_ivec_set (fo,-1,n+1); 636ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 637827bd09bSSatish Balay 6388cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 6398cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 640827bd09bSSatish Balay 641827bd09bSSatish Balay /* determine the # of unique dof */ 642ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 643ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 644ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 645ca8e9878SJed Brown ierr = PetscInfo(0,"done first PCTFS_gs_gop_hc\n");CHKERRQ(ierr); 646ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 647947b95d8SBarry Smith for (i=0; i<n; i++) { 648db4deed7SKarl Rupp if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; } 6492fa5cd67SKarl Rupp if (lhs[i]!=1.0) shared=1; 650827bd09bSSatish Balay } 651827bd09bSSatish Balay 652b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 653827bd09bSSatish Balay rsum[0]+=0.1; 654827bd09bSSatish Balay rsum[1]+=0.1; 655827bd09bSSatish Balay 65652f87cdaSBarry Smith xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 65752f87cdaSBarry Smith xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 658827bd09bSSatish Balay 659827bd09bSSatish Balay /* determine separator sets top down */ 6602fa5cd67SKarl Rupp if (shared) { 661827bd09bSSatish Balay /* solution is to do as in the symmetric shared case but then */ 662827bd09bSSatish Balay /* pick the sub-hc with the most free dofs and do a mat-vec */ 663827bd09bSSatish Balay /* and pick up the responses on the other sub-hc from the */ 664827bd09bSSatish Balay /* initial separator set obtained from the symm. shared case */ 665e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!\n"); 6664e619c20SJed Brown /* [dead code deleted since it is unlikely to be completed] */ 667db4deed7SKarl Rupp } else { 668db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 669827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 670ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 671ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 672827bd09bSSatish Balay 673827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 674ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 675ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 676827bd09bSSatish Balay 677827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 678db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 679db4deed7SKarl Rupp if (!used[i]) { 680827bd09bSSatish Balay /* number of unmarked dofs on node */ 681827bd09bSSatish Balay ct++; 682827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6832fa5cd67SKarl Rupp if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++; 684827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6852fa5cd67SKarl Rupp if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++; 686827bd09bSSatish Balay } 687827bd09bSSatish Balay } 688827bd09bSSatish Balay 689827bd09bSSatish Balay /* for the non-symmetric case we need separators of width 2 */ 690827bd09bSSatish Balay /* so take both sides */ 691827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 692b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 693827bd09bSSatish Balay 694827bd09bSSatish Balay ct=0; 695db4deed7SKarl Rupp if (id<mask) { 696827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 697db4deed7SKarl Rupp for (i=0;i<n;i++) { 698db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 699827bd09bSSatish Balay ct++; nfo++; 700827bd09bSSatish Balay *--iptr = local2global[i]; 701827bd09bSSatish Balay used[i] =edge; 702827bd09bSSatish Balay } 703827bd09bSSatish Balay } 704827bd09bSSatish Balay /* LSH hc summation of ct should be sum[0] */ 705db4deed7SKarl Rupp } else { 706827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 707db4deed7SKarl Rupp for (i=0; i<n; i++) { 708db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 709827bd09bSSatish Balay ct++; nfo++; 710827bd09bSSatish Balay *--iptr = local2global[i]; 711827bd09bSSatish Balay used[i] = edge; 712827bd09bSSatish Balay } 713827bd09bSSatish Balay } 714827bd09bSSatish Balay /* RSH hc summation of ct should be sum[1] */ 715827bd09bSSatish Balay } 716827bd09bSSatish Balay 7172fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 718827bd09bSSatish Balay lnsep[edge]=ct; 719827bd09bSSatish Balay nsep[edge] =sum[0]+sum[1]; 720827bd09bSSatish Balay dir [edge] =BOTH; 721827bd09bSSatish Balay 722827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 723827bd09bSSatish Balay /* do we need full set of separators for this? */ 724827bd09bSSatish Balay 725827bd09bSSatish Balay /* fold rhs hc into lower */ 7262fa5cd67SKarl Rupp if (id>=mask) id-=mask; 727827bd09bSSatish Balay } 728827bd09bSSatish Balay } 729827bd09bSSatish Balay 730827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 731db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 732db4deed7SKarl Rupp if (!used[i]) { 733827bd09bSSatish Balay ct++; nfo++; 734827bd09bSSatish Balay *--iptr = local2global[i]; 735827bd09bSSatish Balay used[i] = edge; 736827bd09bSSatish Balay } 737827bd09bSSatish Balay } 7382fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 739827bd09bSSatish Balay lnsep[edge]=ct; 740827bd09bSSatish Balay nsep [edge]=ct; 741827bd09bSSatish Balay dir [edge]=BOTH; 742827bd09bSSatish Balay 743827bd09bSSatish Balay xyt_handle->info->nsep = nsep; 744827bd09bSSatish Balay xyt_handle->info->lnsep = lnsep; 745827bd09bSSatish Balay xyt_handle->info->fo = fo; 746827bd09bSSatish Balay xyt_handle->info->nfo = nfo; 747827bd09bSSatish Balay 748a501084fSBarry Smith free(dir); 749a501084fSBarry Smith free(lhs); 750a501084fSBarry Smith free(rhs); 751a501084fSBarry Smith free(used); 7520924e98cSBarry Smith PetscFunctionReturn(0); 753827bd09bSSatish Balay } 754827bd09bSSatish Balay 7557b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7565c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data) 757827bd09bSSatish Balay { 758827bd09bSSatish Balay mv_info *mvi; 759827bd09bSSatish Balay 760a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 761827bd09bSSatish Balay mvi->n = n; 762827bd09bSSatish Balay mvi->m = m; 763827bd09bSSatish Balay mvi->n_global = -1; 764827bd09bSSatish Balay mvi->m_global = -1; 76552f87cdaSBarry Smith mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt)); 7662fa5cd67SKarl Rupp 767ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 768827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7695c8f6a95SKarl Rupp mvi->matvec = matvec; 770827bd09bSSatish Balay mvi->grid_data = grid_data; 771827bd09bSSatish Balay 772827bd09bSSatish Balay /* set xyt communication handle to perform restricted matvec */ 773ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 774827bd09bSSatish Balay 775827bd09bSSatish Balay return(mvi); 776827bd09bSSatish Balay } 777827bd09bSSatish Balay 7787b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7793fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 780827bd09bSSatish Balay { 7810924e98cSBarry Smith PetscFunctionBegin; 782827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 7830924e98cSBarry Smith PetscFunctionReturn(0); 784827bd09bSSatish Balay } 785827bd09bSSatish Balay 786827bd09bSSatish Balay 787827bd09bSSatish Balay 788