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 typedef struct matvec_info { 3852f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 3952f87cdaSBarry Smith PetscInt *local2global; 40ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 41a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*); 42827bd09bSSatish Balay void *grid_data; 43827bd09bSSatish Balay } mv_info; 44827bd09bSSatish Balay 45827bd09bSSatish Balay struct xyt_CDT { 4652f87cdaSBarry Smith PetscInt id; 4752f87cdaSBarry Smith PetscInt ns; 4852f87cdaSBarry Smith PetscInt level; 49827bd09bSSatish Balay xyt_info *info; 50827bd09bSSatish Balay mv_info *mvi; 51827bd09bSSatish Balay }; 52827bd09bSSatish Balay 5352f87cdaSBarry Smith static PetscInt n_xyt =0; 5452f87cdaSBarry Smith static PetscInt n_xyt_handles=0; 55827bd09bSSatish Balay 56827bd09bSSatish Balay /* prototypes */ 573fdc5746SBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *rhs); 583fdc5746SBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle); 593fdc5746SBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle); 603fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 6138dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle); 6238dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle); 635c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data); 64827bd09bSSatish Balay 657b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 667b1ae94cSBarry Smith xyt_ADT XYT_new(void) 67827bd09bSSatish Balay { 68827bd09bSSatish Balay xyt_ADT xyt_handle; 69827bd09bSSatish Balay 70827bd09bSSatish Balay /* rolling count on n_xyt ... pot. problem here */ 71827bd09bSSatish Balay n_xyt_handles++; 72a501084fSBarry Smith xyt_handle = (xyt_ADT)malloc(sizeof(struct xyt_CDT)); 73827bd09bSSatish Balay xyt_handle->id = ++n_xyt; 74827bd09bSSatish Balay xyt_handle->info = NULL; 75827bd09bSSatish Balay xyt_handle->mvi = NULL; 76827bd09bSSatish Balay 77827bd09bSSatish Balay return(xyt_handle); 78827bd09bSSatish Balay } 79827bd09bSSatish Balay 807b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 8138dcbb09SBarry Smith PetscErrorCode XYT_factor(xyt_ADT xyt_handle, /* prev. allocated xyt handle */ 8252f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 8352f87cdaSBarry Smith PetscInt n, /* local num rows */ 8452f87cdaSBarry Smith PetscInt m, /* local num cols */ 855c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */ 861147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 87827bd09bSSatish Balay { 88827bd09bSSatish Balay 89b1c944f5SJed Brown PCTFS_comm_init(); 90827bd09bSSatish Balay check_handle(xyt_handle); 91827bd09bSSatish Balay 92827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 932c71b3e2SJacob Faibussowitsch PetscCheckFalse((1<<(xyt_handle->level=PCTFS_i_log2_num_nodes))!=PCTFS_num_nodes,PETSC_COMM_SELF,PETSC_ERR_PLIB,"only 2^k for now and MPI_COMM_WORLD!!! %D != %D",1<<PCTFS_i_log2_num_nodes,PCTFS_num_nodes); 94827bd09bSSatish Balay 95827bd09bSSatish Balay /* space for X info */ 96a501084fSBarry Smith xyt_handle->info = (xyt_info*)malloc(sizeof(xyt_info)); 97827bd09bSSatish Balay 98827bd09bSSatish Balay /* set up matvec handles */ 995c8f6a95SKarl Rupp xyt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data); 100827bd09bSSatish Balay 101827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 102827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 103827bd09bSSatish Balay xyt_handle->ns=0; 104827bd09bSSatish Balay 105827bd09bSSatish Balay /* determine separators and generate firing order - NB xyt info set here */ 106827bd09bSSatish Balay det_separators(xyt_handle); 107827bd09bSSatish Balay 108827bd09bSSatish Balay return(do_xyt_factor(xyt_handle)); 109827bd09bSSatish Balay } 110827bd09bSSatish Balay 1117b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 11238dcbb09SBarry Smith PetscErrorCode XYT_solve(xyt_ADT xyt_handle, PetscScalar *x, PetscScalar *b) 113827bd09bSSatish Balay { 114b1c944f5SJed Brown PCTFS_comm_init(); 115827bd09bSSatish Balay check_handle(xyt_handle); 116827bd09bSSatish Balay 117827bd09bSSatish Balay /* need to copy b into x? */ 11859bc5b24SSatish Balay if (b) PCTFS_rvec_copy(x,b,xyt_handle->mvi->n); 11938dcbb09SBarry Smith return do_xyt_solve(xyt_handle,x); 120827bd09bSSatish Balay } 121827bd09bSSatish Balay 1227b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 12338dcbb09SBarry Smith PetscErrorCode XYT_free(xyt_ADT xyt_handle) 124827bd09bSSatish Balay { 125b1c944f5SJed Brown PCTFS_comm_init(); 126827bd09bSSatish Balay check_handle(xyt_handle); 127827bd09bSSatish Balay n_xyt_handles--; 128827bd09bSSatish Balay 129a501084fSBarry Smith free(xyt_handle->info->nsep); 130a501084fSBarry Smith free(xyt_handle->info->lnsep); 131a501084fSBarry Smith free(xyt_handle->info->fo); 132a501084fSBarry Smith free(xyt_handle->info->stages); 133a501084fSBarry Smith free(xyt_handle->info->solve_uu); 134a501084fSBarry Smith free(xyt_handle->info->solve_w); 135a501084fSBarry Smith free(xyt_handle->info->x); 136a501084fSBarry Smith free(xyt_handle->info->xcol_vals); 137a501084fSBarry Smith free(xyt_handle->info->xcol_sz); 138a501084fSBarry Smith free(xyt_handle->info->xcol_indices); 139a501084fSBarry Smith free(xyt_handle->info->y); 140a501084fSBarry Smith free(xyt_handle->info->ycol_vals); 141a501084fSBarry Smith free(xyt_handle->info->ycol_sz); 142a501084fSBarry Smith free(xyt_handle->info->ycol_indices); 143a501084fSBarry Smith free(xyt_handle->info); 144a501084fSBarry Smith free(xyt_handle->mvi->local2global); 145ca8e9878SJed Brown PCTFS_gs_free(xyt_handle->mvi->PCTFS_gs_handle); 146a501084fSBarry Smith free(xyt_handle->mvi); 147a501084fSBarry Smith free(xyt_handle); 148827bd09bSSatish Balay 149827bd09bSSatish Balay /* if the check fails we nuke */ 150a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 151827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 152827bd09bSSatish Balay return(0); 153827bd09bSSatish Balay } 154827bd09bSSatish Balay 1557b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 15638dcbb09SBarry Smith PetscErrorCode XYT_stats(xyt_ADT xyt_handle) 157827bd09bSSatish Balay { 15852f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 15952f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 16052f87cdaSBarry Smith PetscInt vals[9], work[9]; 161a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 162827bd09bSSatish Balay 163b1c944f5SJed Brown PCTFS_comm_init(); 164827bd09bSSatish Balay check_handle(xyt_handle); 165827bd09bSSatish Balay 166827bd09bSSatish Balay /* if factorization not done there are no stats */ 1677d0a6c19SBarry Smith if (!xyt_handle->info||!xyt_handle->mvi) { 168b1c944f5SJed Brown if (!PCTFS_my_id) PetscPrintf(PETSC_COMM_WORLD,"XYT_stats() :: no stats available!\n"); 169827bd09bSSatish Balay return 1; 170827bd09bSSatish Balay } 171827bd09bSSatish Balay 172827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xyt_handle->info->nnz; 173827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xyt_handle->mvi->n; 174827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xyt_handle->info->msg_buf_sz; 175*dd39110bSPierre Jolivet PCTFS_giop(vals,work,PETSC_STATIC_ARRAY_LENGTH(op)-1,op); 176827bd09bSSatish Balay 1772fa5cd67SKarl Rupp fvals[0]=fvals[1]=fvals[2]=xyt_handle->info->tot_solve_time/xyt_handle->info->nsolves++; 178*dd39110bSPierre Jolivet PCTFS_grop(fvals,fwork,PETSC_STATIC_ARRAY_LENGTH(fop)-1,fop); 179827bd09bSSatish Balay 180b1c944f5SJed Brown if (!PCTFS_my_id) { 181b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_nnz=%D\n",PCTFS_my_id,vals[0]); 182b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_nnz=%D\n",PCTFS_my_id,vals[1]); 183b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes); 184b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_nnz=%D\n",PCTFS_my_id,vals[2]); 18577b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(2d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5))); 18677b4d14cSPeter Brune PetscPrintf(PETSC_COMM_WORLD,"%D :: xyt C(3d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667))); 187b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_n =%D\n",PCTFS_my_id,vals[3]); 188b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_n =%D\n",PCTFS_my_id,vals[4]); 189b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes); 190b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xyt_n =%D\n",PCTFS_my_id,vals[5]); 191b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_buf=%D\n",PCTFS_my_id,vals[6]); 192b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_buf=%D\n",PCTFS_my_id,vals[7]); 193b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes); 194b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: min xyt_slv=%g\n",PCTFS_my_id,fvals[0]); 195b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: max xyt_slv=%g\n",PCTFS_my_id,fvals[1]); 196b1c944f5SJed Brown PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xyt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes); 197827bd09bSSatish Balay } 198827bd09bSSatish Balay 199827bd09bSSatish Balay return(0); 200827bd09bSSatish Balay } 201827bd09bSSatish Balay 202827bd09bSSatish Balay /*************************************xyt.c************************************ 203827bd09bSSatish Balay 204827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 205827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 206827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 207827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 208827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 209827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 210ca8e9878SJed Brown PCTFS_gs_init/gop). 211827bd09bSSatish Balay 212827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 213827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 214827bd09bSSatish Balay **************************************xyt.c***********************************/ 21538dcbb09SBarry Smith static PetscErrorCode do_xyt_factor(xyt_ADT xyt_handle) 216827bd09bSSatish Balay { 2177b1ae94cSBarry Smith return xyt_generate(xyt_handle); 218827bd09bSSatish Balay } 219827bd09bSSatish Balay 2207b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 22138dcbb09SBarry Smith static PetscErrorCode xyt_generate(xyt_ADT xyt_handle) 222827bd09bSSatish Balay { 22352f87cdaSBarry Smith PetscInt i,j,k,idx; 22452f87cdaSBarry Smith PetscInt dim, col; 225a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22652f87cdaSBarry Smith PetscInt *segs; 22752f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 22852f87cdaSBarry Smith PetscInt off, len; 229a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr; 23052f87cdaSBarry Smith PetscInt *iptr, flag; 23152f87cdaSBarry Smith PetscInt start =0, end, work; 23252f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 233ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 23452f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 23552f87cdaSBarry Smith PetscInt a_n =xyt_handle->mvi->n; 23652f87cdaSBarry Smith PetscInt a_m =xyt_handle->mvi->m; 23752f87cdaSBarry Smith PetscInt *a_local2global=xyt_handle->mvi->local2global; 23852f87cdaSBarry Smith PetscInt level; 23952f87cdaSBarry Smith PetscInt n, m; 24052f87cdaSBarry Smith PetscInt *xcol_sz, *xcol_indices, *stages; 241a501084fSBarry Smith PetscScalar **xcol_vals, *x; 24252f87cdaSBarry Smith PetscInt *ycol_sz, *ycol_indices; 243a501084fSBarry Smith PetscScalar **ycol_vals, *y; 24452f87cdaSBarry Smith PetscInt n_global; 24552f87cdaSBarry Smith PetscInt xt_nnz =0, xt_max_nnz=0; 24652f87cdaSBarry Smith PetscInt yt_nnz =0, yt_max_nnz=0; 24752f87cdaSBarry Smith PetscInt xt_zero_nnz =0; 24852f87cdaSBarry Smith PetscInt xt_zero_nnz_0=0; 24952f87cdaSBarry Smith PetscInt yt_zero_nnz =0; 25052f87cdaSBarry Smith PetscInt yt_zero_nnz_0=0; 2514a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 252a501084fSBarry Smith PetscScalar dm1 = -1.0; 253827bd09bSSatish Balay 254827bd09bSSatish Balay n =xyt_handle->mvi->n; 255827bd09bSSatish Balay nsep =xyt_handle->info->nsep; 256827bd09bSSatish Balay lnsep =xyt_handle->info->lnsep; 257827bd09bSSatish Balay fo =xyt_handle->info->fo; 258827bd09bSSatish Balay end =lnsep[0]; 259827bd09bSSatish Balay level =xyt_handle->level; 260ca8e9878SJed Brown PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 261827bd09bSSatish Balay 262827bd09bSSatish Balay /* is there a null space? */ 263827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2642fa5cd67SKarl Rupp for (i=0, j=0; i<=level; i++) j+=nsep[i]; 265827bd09bSSatish Balay 266827bd09bSSatish Balay m = j-xyt_handle->ns; 26722d28d08SBarry Smith if (m!=j) { 2689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"xyt_generate() :: null space exists %D %D %D\n",m,j,xyt_handle->ns)); 26922d28d08SBarry Smith } 270827bd09bSSatish Balay 2719566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"xyt_generate() :: X(%D,%D)\n",n,m)); 272827bd09bSSatish Balay 273827bd09bSSatish Balay /* get and initialize storage for x local */ 274827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 27552f87cdaSBarry Smith xcol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 27652f87cdaSBarry Smith xcol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 277a501084fSBarry Smith xcol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 278db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 279827bd09bSSatish Balay xcol_indices[j]=xcol_indices[j+1]=xcol_sz[i]=-1; 280827bd09bSSatish Balay xcol_vals[i] = NULL; 281827bd09bSSatish Balay } 282827bd09bSSatish Balay xcol_indices[j]=-1; 283827bd09bSSatish Balay 284827bd09bSSatish Balay /* get and initialize storage for y local */ 285827bd09bSSatish Balay /* note that y local is nxm and stored by columns */ 28652f87cdaSBarry Smith ycol_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 28752f87cdaSBarry Smith ycol_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 288a501084fSBarry Smith ycol_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 289db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 290827bd09bSSatish Balay ycol_indices[j]=ycol_indices[j+1]=ycol_sz[i]=-1; 291827bd09bSSatish Balay ycol_vals[i] = NULL; 292827bd09bSSatish Balay } 293827bd09bSSatish Balay ycol_indices[j]=-1; 294827bd09bSSatish Balay 295827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 296827bd09bSSatish Balay /* this looks like nsep[]=segments */ 29752f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 29852f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 299ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 300ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 3012fa5cd67SKarl Rupp for (i=0; i<level; i++) segs[i+1] += segs[i]; 302827bd09bSSatish Balay stages[0] = segs[0]; 303827bd09bSSatish Balay 304827bd09bSSatish Balay /* temporary vectors */ 305a501084fSBarry Smith u = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 306a501084fSBarry Smith z = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 307a501084fSBarry Smith v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar)); 308a501084fSBarry Smith uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 309a501084fSBarry Smith w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 310827bd09bSSatish Balay 311827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 3122fa5cd67SKarl Rupp for (i=1, j=0; i<=level; i++) j+=nsep[i]; 313827bd09bSSatish Balay 314827bd09bSSatish Balay /* storage for sparse x values */ 315827bd09bSSatish Balay n_global = xyt_handle->info->n_global; 31685ec1a3cSBarry Smith xt_max_nnz = yt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 317a501084fSBarry Smith x = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 318a501084fSBarry Smith y = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 319827bd09bSSatish Balay 320827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 321827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 322db4deed7SKarl Rupp for (dim=i=j=0; i<m; i++) { 323827bd09bSSatish Balay /* time to move to the next level? */ 324d4af0d30SBarry Smith while (i==segs[dim]) { 3252c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim==level,PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level"); 326827bd09bSSatish Balay stages[dim++]=i; 327827bd09bSSatish Balay end +=lnsep[dim]; 328827bd09bSSatish Balay } 329827bd09bSSatish Balay stages[dim]=i; 330827bd09bSSatish Balay 331827bd09bSSatish Balay /* which column are we firing? */ 332827bd09bSSatish Balay /* i.e. set v_l */ 333827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 334827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 335b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 336827bd09bSSatish Balay 337827bd09bSSatish Balay /* shouldn't need this */ 338db4deed7SKarl Rupp if (col==INT_MAX) { 3399566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"hey ... col==INT_MAX??\n")); 340827bd09bSSatish Balay continue; 341827bd09bSSatish Balay } 342827bd09bSSatish Balay 343827bd09bSSatish Balay /* do I own it? I should */ 344ca8e9878SJed Brown PCTFS_rvec_zero(v,a_m); 3452fa5cd67SKarl Rupp if (col==fo[start]) { 346827bd09bSSatish Balay start++; 347ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_n); 3482fa5cd67SKarl Rupp if (idx!=-1) { 3492fa5cd67SKarl Rupp v[idx] = 1.0; 3502fa5cd67SKarl Rupp j++; 351546078acSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!"); 352db4deed7SKarl Rupp } else { 353ca8e9878SJed Brown idx=PCTFS_ivec_linear_search(col, a_local2global, a_m); 3542fa5cd67SKarl Rupp if (idx!=-1) v[idx] = 1.0; 355827bd09bSSatish Balay } 356827bd09bSSatish Balay 357827bd09bSSatish Balay /* perform u = A.v_l */ 358ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 359827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 360827bd09bSSatish Balay 361827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 362827bd09bSSatish Balay /* technically only need to zero out first i entries */ 363827bd09bSSatish Balay /* later turn this into an XYT_solve call ? */ 364ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 365827bd09bSSatish Balay y_ptr=y; 366827bd09bSSatish Balay iptr = ycol_indices; 367db4deed7SKarl Rupp for (k=0; k<i; k++) { 368827bd09bSSatish Balay off = *iptr++; 369827bd09bSSatish Balay len = *iptr++; 3709566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len,&dlen)); 3718b83055fSJed Brown PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,y_ptr,&i1)); 372827bd09bSSatish Balay y_ptr+=len; 373827bd09bSSatish Balay } 374827bd09bSSatish Balay 375827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 3769566063dSJacob Faibussowitsch PetscCall(PCTFS_ssgl_radd (uu, w, dim, stages)); 377827bd09bSSatish Balay 378827bd09bSSatish Balay /* z = X.uu */ 379ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 380827bd09bSSatish Balay x_ptr=x; 381827bd09bSSatish Balay iptr = xcol_indices; 382db4deed7SKarl Rupp for (k=0; k<i; k++) { 383827bd09bSSatish Balay off = *iptr++; 384827bd09bSSatish Balay len = *iptr++; 3859566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len,&dlen)); 3868b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1)); 387827bd09bSSatish Balay x_ptr+=len; 388827bd09bSSatish Balay } 389827bd09bSSatish Balay 390827bd09bSSatish Balay /* compute v_l = v_l - z */ 391ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 3929566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n,&dlen)); 3938b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1)); 394827bd09bSSatish Balay 395827bd09bSSatish Balay /* compute u_l = A.v_l */ 3962fa5cd67SKarl Rupp if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); 397ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 398827bd09bSSatish Balay do_matvec(xyt_handle->mvi,v,u); 399827bd09bSSatish Balay 400827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - local portion */ 4019566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n,&dlen)); 4028b83055fSJed Brown PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,u,&i1)); 403827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(u_l^T.u_l) - comm portion */ 404b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 405827bd09bSSatish Balay 4068f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 407827bd09bSSatish Balay 408827bd09bSSatish Balay /* check for small alpha */ 409827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 4102c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsScalar(alpha)<1.0e-14,PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g",alpha); 411827bd09bSSatish Balay 412827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 413ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 414ca8e9878SJed Brown PCTFS_rvec_scale(u,1.0/alpha,n); 415827bd09bSSatish Balay 416827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 417827bd09bSSatish Balay flag = 1; 418827bd09bSSatish Balay off =len=0; 419db4deed7SKarl Rupp for (k=0; k<n; k++) { 420db4deed7SKarl Rupp if (v[k]!=0.0) { 421827bd09bSSatish Balay len=k; 422db4deed7SKarl Rupp if (flag) {off=k; flag=0;} 423827bd09bSSatish Balay } 424827bd09bSSatish Balay } 425827bd09bSSatish Balay 426827bd09bSSatish Balay len -= (off-1); 427827bd09bSSatish Balay 428db4deed7SKarl Rupp if (len>0) { 429db4deed7SKarl Rupp if ((xt_nnz+len)>xt_max_nnz) { 4309566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"increasing space for X by 2x!\n")); 431827bd09bSSatish Balay xt_max_nnz *= 2; 432a501084fSBarry Smith x_ptr = (PetscScalar*) malloc(xt_max_nnz*sizeof(PetscScalar)); 433ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xt_nnz); 434a501084fSBarry Smith free(x); 435827bd09bSSatish Balay x = x_ptr; 436827bd09bSSatish Balay x_ptr+=xt_nnz; 437827bd09bSSatish Balay } 438827bd09bSSatish Balay xt_nnz += len; 439ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 440827bd09bSSatish Balay 441827bd09bSSatish Balay /* keep track of number of zeros */ 442db4deed7SKarl Rupp if (dim) { 443db4deed7SKarl Rupp for (k=0; k<len; k++) { 4442fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz++; 445827bd09bSSatish Balay } 446db4deed7SKarl Rupp } else { 447db4deed7SKarl Rupp for (k=0; k<len; k++) { 4482fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xt_zero_nnz_0++; 449827bd09bSSatish Balay } 450827bd09bSSatish Balay } 451827bd09bSSatish Balay xcol_indices[2*i] = off; 452827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = len; 453827bd09bSSatish Balay xcol_vals[i] = x_ptr; 454db4deed7SKarl Rupp } else { 455827bd09bSSatish Balay xcol_indices[2*i] = 0; 456827bd09bSSatish Balay xcol_sz[i] = xcol_indices[2*i+1] = 0; 457827bd09bSSatish Balay xcol_vals[i] = x_ptr; 458827bd09bSSatish Balay } 459827bd09bSSatish Balay 460827bd09bSSatish Balay /* add newly generated column, u_l, to Y */ 461827bd09bSSatish Balay flag = 1; 462827bd09bSSatish Balay off =len=0; 463db4deed7SKarl Rupp for (k=0; k<n; k++) { 464db4deed7SKarl Rupp if (u[k]!=0.0) { 465827bd09bSSatish Balay len=k; 466db4deed7SKarl Rupp if (flag) { off=k; flag=0; } 467827bd09bSSatish Balay } 468827bd09bSSatish Balay } 469827bd09bSSatish Balay 470827bd09bSSatish Balay len -= (off-1); 471827bd09bSSatish Balay 472db4deed7SKarl Rupp if (len>0) { 473db4deed7SKarl Rupp if ((yt_nnz+len)>yt_max_nnz) { 4749566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"increasing space for Y by 2x!\n")); 475827bd09bSSatish Balay yt_max_nnz *= 2; 476a501084fSBarry Smith y_ptr = (PetscScalar*) malloc(yt_max_nnz*sizeof(PetscScalar)); 477ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,y,yt_nnz); 478a501084fSBarry Smith free(y); 479827bd09bSSatish Balay y = y_ptr; 480827bd09bSSatish Balay y_ptr+=yt_nnz; 481827bd09bSSatish Balay } 482827bd09bSSatish Balay yt_nnz += len; 483ca8e9878SJed Brown PCTFS_rvec_copy(y_ptr,u+off,len); 484827bd09bSSatish Balay 485827bd09bSSatish Balay /* keep track of number of zeros */ 486db4deed7SKarl Rupp if (dim) { 487db4deed7SKarl Rupp for (k=0; k<len; k++) { 4882fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz++; 489827bd09bSSatish Balay } 490db4deed7SKarl Rupp } else { 491db4deed7SKarl Rupp for (k=0; k<len; k++) { 4922fa5cd67SKarl Rupp if (y_ptr[k]==0.0) yt_zero_nnz_0++; 493827bd09bSSatish Balay } 494827bd09bSSatish Balay } 495827bd09bSSatish Balay ycol_indices[2*i] = off; 496827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = len; 497827bd09bSSatish Balay ycol_vals[i] = y_ptr; 498db4deed7SKarl Rupp } else { 499827bd09bSSatish Balay ycol_indices[2*i] = 0; 500827bd09bSSatish Balay ycol_sz[i] = ycol_indices[2*i+1] = 0; 501827bd09bSSatish Balay ycol_vals[i] = y_ptr; 502827bd09bSSatish Balay } 503827bd09bSSatish Balay } 504827bd09bSSatish Balay 505827bd09bSSatish Balay /* close off stages for execution phase */ 506db4deed7SKarl Rupp while (dim!=level) { 507827bd09bSSatish Balay stages[dim++]=i; 5089566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level)); 509827bd09bSSatish Balay } 510827bd09bSSatish Balay stages[dim]=i; 511827bd09bSSatish Balay 512827bd09bSSatish Balay xyt_handle->info->n =xyt_handle->mvi->n; 513827bd09bSSatish Balay xyt_handle->info->m =m; 514827bd09bSSatish Balay xyt_handle->info->nnz =xt_nnz + yt_nnz; 515827bd09bSSatish Balay xyt_handle->info->max_nnz =xt_max_nnz + yt_max_nnz; 516827bd09bSSatish Balay xyt_handle->info->msg_buf_sz =stages[level]-stages[0]; 517a501084fSBarry Smith xyt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 518a501084fSBarry Smith xyt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 519827bd09bSSatish Balay xyt_handle->info->x =x; 520827bd09bSSatish Balay xyt_handle->info->xcol_vals =xcol_vals; 521827bd09bSSatish Balay xyt_handle->info->xcol_sz =xcol_sz; 522827bd09bSSatish Balay xyt_handle->info->xcol_indices=xcol_indices; 523827bd09bSSatish Balay xyt_handle->info->stages =stages; 524827bd09bSSatish Balay xyt_handle->info->y =y; 525827bd09bSSatish Balay xyt_handle->info->ycol_vals =ycol_vals; 526827bd09bSSatish Balay xyt_handle->info->ycol_sz =ycol_sz; 527827bd09bSSatish Balay xyt_handle->info->ycol_indices=ycol_indices; 528827bd09bSSatish Balay 529a501084fSBarry Smith free(segs); 530a501084fSBarry Smith free(u); 531a501084fSBarry Smith free(v); 532a501084fSBarry Smith free(uu); 533a501084fSBarry Smith free(z); 534a501084fSBarry Smith free(w); 535827bd09bSSatish Balay 536827bd09bSSatish Balay return(0); 537827bd09bSSatish Balay } 538827bd09bSSatish Balay 5397b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5400924e98cSBarry Smith static PetscErrorCode do_xyt_solve(xyt_ADT xyt_handle, PetscScalar *uc) 541827bd09bSSatish Balay { 54252f87cdaSBarry Smith PetscInt off, len, *iptr; 54352f87cdaSBarry Smith PetscInt level =xyt_handle->level; 54452f87cdaSBarry Smith PetscInt n =xyt_handle->info->n; 54552f87cdaSBarry Smith PetscInt m =xyt_handle->info->m; 54652f87cdaSBarry Smith PetscInt *stages =xyt_handle->info->stages; 54752f87cdaSBarry Smith PetscInt *xcol_indices=xyt_handle->info->xcol_indices; 54852f87cdaSBarry Smith PetscInt *ycol_indices=xyt_handle->info->ycol_indices; 549a501084fSBarry Smith PetscScalar *x_ptr, *y_ptr, *uu_ptr; 550a501084fSBarry Smith PetscScalar *solve_uu=xyt_handle->info->solve_uu; 551a501084fSBarry Smith PetscScalar *solve_w =xyt_handle->info->solve_w; 552a501084fSBarry Smith PetscScalar *x =xyt_handle->info->x; 553a501084fSBarry Smith PetscScalar *y =xyt_handle->info->y; 5544a0de3f6SBarry Smith PetscBLASInt i1 = 1,dlen; 555827bd09bSSatish Balay 5560924e98cSBarry Smith PetscFunctionBegin; 557827bd09bSSatish Balay uu_ptr=solve_uu; 558ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 559827bd09bSSatish Balay 560827bd09bSSatish Balay /* x = X.Y^T.b */ 561827bd09bSSatish Balay /* uu = Y^T.b */ 562947b95d8SBarry Smith for (y_ptr=y,iptr=ycol_indices; *iptr!=-1; y_ptr+=len) { 563c5df96a5SBarry Smith off =*iptr++; 564c5df96a5SBarry Smith len =*iptr++; 5659566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len,&dlen)); 5668b83055fSJed Brown PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,y_ptr,&i1)); 567827bd09bSSatish Balay } 568827bd09bSSatish Balay 569827bd09bSSatish Balay /* comunication of beta */ 570827bd09bSSatish Balay uu_ptr=solve_uu; 5719566063dSJacob Faibussowitsch if (level) PetscCall(PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages)); 572ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 573827bd09bSSatish Balay 574827bd09bSSatish Balay /* x = X.uu */ 575db4deed7SKarl Rupp for (x_ptr=x,iptr=xcol_indices; *iptr!=-1; x_ptr+=len) { 576c5df96a5SBarry Smith off =*iptr++; 577c5df96a5SBarry Smith len =*iptr++; 5789566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(len,&dlen)); 5798b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1)); 580827bd09bSSatish Balay } 5810924e98cSBarry Smith PetscFunctionReturn(0); 582827bd09bSSatish Balay } 583827bd09bSSatish Balay 5847b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5850924e98cSBarry Smith static PetscErrorCode check_handle(xyt_ADT xyt_handle) 586827bd09bSSatish Balay { 58752f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 588827bd09bSSatish Balay 5890924e98cSBarry Smith PetscFunctionBegin; 59028b400f6SJacob Faibussowitsch PetscCheck(xyt_handle,PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D",xyt_handle); 591827bd09bSSatish Balay 592827bd09bSSatish Balay vals[0]=vals[1]=xyt_handle->id; 593*dd39110bSPierre Jolivet PCTFS_giop(vals,work,PETSC_STATIC_ARRAY_LENGTH(op)-1,op); 5942c71b3e2SJacob Faibussowitsch PetscCheckFalse((vals[0]!=vals[1])||(xyt_handle->id<=0),PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: id mismatch min/max %D/%D %D", vals[0],vals[1], xyt_handle->id); 5950924e98cSBarry Smith PetscFunctionReturn(0); 596827bd09bSSatish Balay } 597827bd09bSSatish Balay 5987b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 5990924e98cSBarry Smith static PetscErrorCode det_separators(xyt_ADT xyt_handle) 600827bd09bSSatish Balay { 60152f87cdaSBarry Smith PetscInt i, ct, id; 60252f87cdaSBarry Smith PetscInt mask, edge, *iptr; 60352f87cdaSBarry Smith PetscInt *dir, *used; 60452f87cdaSBarry Smith PetscInt sum[4], w[4]; 605a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 60652f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 607a501084fSBarry Smith PetscScalar *lhs, *rhs; 60852f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 609ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle=xyt_handle->mvi->PCTFS_gs_handle; 61052f87cdaSBarry Smith PetscInt *local2global =xyt_handle->mvi->local2global; 61152f87cdaSBarry Smith PetscInt n =xyt_handle->mvi->n; 61252f87cdaSBarry Smith PetscInt m =xyt_handle->mvi->m; 61352f87cdaSBarry Smith PetscInt level =xyt_handle->level; 614ab824b78SBarry Smith PetscInt shared =0; 615827bd09bSSatish Balay 6160924e98cSBarry Smith PetscFunctionBegin; 61752f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 61852f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 61952f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 62052f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 62152f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 622827bd09bSSatish Balay 623ca8e9878SJed Brown PCTFS_ivec_zero(dir,level+1); 624ca8e9878SJed Brown PCTFS_ivec_zero(nsep,level+1); 625ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 626ca8e9878SJed Brown PCTFS_ivec_set (fo,-1,n+1); 627ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 628827bd09bSSatish Balay 6298cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 6308cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 631827bd09bSSatish Balay 632827bd09bSSatish Balay /* determine the # of unique dof */ 633ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 634ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 635ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 6369566063dSJacob Faibussowitsch PetscCall(PetscInfo(0,"done first PCTFS_gs_gop_hc\n")); 637ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 638947b95d8SBarry Smith for (i=0; i<n; i++) { 639db4deed7SKarl Rupp if (lhs[i]!=0.0) { rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; } 6402fa5cd67SKarl Rupp if (lhs[i]!=1.0) shared=1; 641827bd09bSSatish Balay } 642827bd09bSSatish Balay 643b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 644827bd09bSSatish Balay rsum[0]+=0.1; 645827bd09bSSatish Balay rsum[1]+=0.1; 646827bd09bSSatish Balay 64752f87cdaSBarry Smith xyt_handle->info->n_global=xyt_handle->info->m_global=(PetscInt) rsum[0]; 64852f87cdaSBarry Smith xyt_handle->mvi->n_global =xyt_handle->mvi->m_global =(PetscInt) rsum[0]; 649827bd09bSSatish Balay 650827bd09bSSatish Balay /* determine separator sets top down */ 6512fa5cd67SKarl Rupp if (shared) { 652827bd09bSSatish Balay /* solution is to do as in the symmetric shared case but then */ 653827bd09bSSatish Balay /* pick the sub-hc with the most free dofs and do a mat-vec */ 654827bd09bSSatish Balay /* and pick up the responses on the other sub-hc from the */ 655827bd09bSSatish Balay /* initial separator set obtained from the symm. shared case */ 656546078acSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"shared dof separator determination not ready ... see hmt!!!"); 6574e619c20SJed Brown /* [dead code deleted since it is unlikely to be completed] */ 658db4deed7SKarl Rupp } else { 659db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 660827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 661ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 662ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 663827bd09bSSatish Balay 664827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 665ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 666ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 667827bd09bSSatish Balay 668827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 669db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 670db4deed7SKarl Rupp if (!used[i]) { 671827bd09bSSatish Balay /* number of unmarked dofs on node */ 672827bd09bSSatish Balay ct++; 673827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6742fa5cd67SKarl Rupp if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++; 675827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6762fa5cd67SKarl Rupp if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++; 677827bd09bSSatish Balay } 678827bd09bSSatish Balay } 679827bd09bSSatish Balay 680827bd09bSSatish Balay /* for the non-symmetric case we need separators of width 2 */ 681827bd09bSSatish Balay /* so take both sides */ 682827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 683b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 684827bd09bSSatish Balay 685827bd09bSSatish Balay ct=0; 686db4deed7SKarl Rupp if (id<mask) { 687827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 688db4deed7SKarl Rupp for (i=0;i<n;i++) { 689db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 690827bd09bSSatish Balay ct++; nfo++; 691827bd09bSSatish Balay *--iptr = local2global[i]; 692827bd09bSSatish Balay used[i] =edge; 693827bd09bSSatish Balay } 694827bd09bSSatish Balay } 695827bd09bSSatish Balay /* LSH hc summation of ct should be sum[0] */ 696db4deed7SKarl Rupp } else { 697827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 698db4deed7SKarl Rupp for (i=0; i<n; i++) { 699db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 700827bd09bSSatish Balay ct++; nfo++; 701827bd09bSSatish Balay *--iptr = local2global[i]; 702827bd09bSSatish Balay used[i] = edge; 703827bd09bSSatish Balay } 704827bd09bSSatish Balay } 705827bd09bSSatish Balay /* RSH hc summation of ct should be sum[1] */ 706827bd09bSSatish Balay } 707827bd09bSSatish Balay 7082fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 709827bd09bSSatish Balay lnsep[edge]=ct; 710827bd09bSSatish Balay nsep[edge] =sum[0]+sum[1]; 711827bd09bSSatish Balay dir [edge] =BOTH; 712827bd09bSSatish Balay 713827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 714827bd09bSSatish Balay /* do we need full set of separators for this? */ 715827bd09bSSatish Balay 716827bd09bSSatish Balay /* fold rhs hc into lower */ 7172fa5cd67SKarl Rupp if (id>=mask) id-=mask; 718827bd09bSSatish Balay } 719827bd09bSSatish Balay } 720827bd09bSSatish Balay 721827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 722db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 723db4deed7SKarl Rupp if (!used[i]) { 724827bd09bSSatish Balay ct++; nfo++; 725827bd09bSSatish Balay *--iptr = local2global[i]; 726827bd09bSSatish Balay used[i] = edge; 727827bd09bSSatish Balay } 728827bd09bSSatish Balay } 7292fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 730827bd09bSSatish Balay lnsep[edge]=ct; 731827bd09bSSatish Balay nsep [edge]=ct; 732827bd09bSSatish Balay dir [edge]=BOTH; 733827bd09bSSatish Balay 734827bd09bSSatish Balay xyt_handle->info->nsep = nsep; 735827bd09bSSatish Balay xyt_handle->info->lnsep = lnsep; 736827bd09bSSatish Balay xyt_handle->info->fo = fo; 737827bd09bSSatish Balay xyt_handle->info->nfo = nfo; 738827bd09bSSatish Balay 739a501084fSBarry Smith free(dir); 740a501084fSBarry Smith free(lhs); 741a501084fSBarry Smith free(rhs); 742a501084fSBarry Smith free(used); 7430924e98cSBarry Smith PetscFunctionReturn(0); 744827bd09bSSatish Balay } 745827bd09bSSatish Balay 7467b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7475c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), void *grid_data) 748827bd09bSSatish Balay { 749827bd09bSSatish Balay mv_info *mvi; 750827bd09bSSatish Balay 751a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 752827bd09bSSatish Balay mvi->n = n; 753827bd09bSSatish Balay mvi->m = m; 754827bd09bSSatish Balay mvi->n_global = -1; 755827bd09bSSatish Balay mvi->m_global = -1; 75652f87cdaSBarry Smith mvi->local2global= (PetscInt*)malloc((m+1)*sizeof(PetscInt)); 7572fa5cd67SKarl Rupp 758ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 759827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7605c8f6a95SKarl Rupp mvi->matvec = matvec; 761827bd09bSSatish Balay mvi->grid_data = grid_data; 762827bd09bSSatish Balay 763827bd09bSSatish Balay /* set xyt communication handle to perform restricted matvec */ 764ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 765827bd09bSSatish Balay 766827bd09bSSatish Balay return(mvi); 767827bd09bSSatish Balay } 768827bd09bSSatish Balay 7697b1ae94cSBarry Smith /**************************************xyt.c***********************************/ 7703fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 771827bd09bSSatish Balay { 7720924e98cSBarry Smith PetscFunctionBegin; 773827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 7740924e98cSBarry Smith PetscFunctionReturn(0); 775827bd09bSSatish Balay } 776