1827bd09bSSatish Balay 2827bd09bSSatish Balay /*************************************xxt.c************************************ 3827bd09bSSatish Balay Module Name: xxt 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 **************************************xxt.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 xxt_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 *col_sz, *col_indices; 30a501084fSBarry Smith PetscScalar **col_vals, *x, *solve_uu, *solve_w; 3152f87cdaSBarry Smith PetscInt nsolves; 32a501084fSBarry Smith PetscScalar tot_solve_time; 33827bd09bSSatish Balay } xxt_info; 34827bd09bSSatish Balay 35827bd09bSSatish Balay typedef struct matvec_info { 3652f87cdaSBarry Smith PetscInt n, m, n_global, m_global; 3752f87cdaSBarry Smith PetscInt *local2global; 38*ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 39a501084fSBarry Smith PetscErrorCode (*matvec)(struct matvec_info*,PetscScalar*,PetscScalar*); 40827bd09bSSatish Balay void *grid_data; 41827bd09bSSatish Balay } mv_info; 42827bd09bSSatish Balay 43827bd09bSSatish Balay struct xxt_CDT{ 4452f87cdaSBarry Smith PetscInt id; 4552f87cdaSBarry Smith PetscInt ns; 4652f87cdaSBarry Smith PetscInt level; 47827bd09bSSatish Balay xxt_info *info; 48827bd09bSSatish Balay mv_info *mvi; 49827bd09bSSatish Balay }; 50827bd09bSSatish Balay 5152f87cdaSBarry Smith static PetscInt n_xxt=0; 5252f87cdaSBarry Smith static PetscInt n_xxt_handles=0; 53827bd09bSSatish Balay 54827bd09bSSatish Balay /* prototypes */ 553fdc5746SBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *rhs); 563fdc5746SBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle); 573fdc5746SBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle); 583fdc5746SBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u); 5952f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle); 6052f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle); 6152f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data); 62a501084fSBarry Smith 637b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 647b1ae94cSBarry Smith xxt_ADT XXT_new(void) 65827bd09bSSatish Balay { 66827bd09bSSatish Balay xxt_ADT xxt_handle; 67827bd09bSSatish Balay 68827bd09bSSatish Balay /* rolling count on n_xxt ... pot. problem here */ 69827bd09bSSatish Balay n_xxt_handles++; 70a501084fSBarry Smith xxt_handle = (xxt_ADT)malloc(sizeof(struct xxt_CDT)); 71827bd09bSSatish Balay xxt_handle->id = ++n_xxt; 72827bd09bSSatish Balay xxt_handle->info = NULL; xxt_handle->mvi = NULL; 73827bd09bSSatish Balay 74827bd09bSSatish Balay return(xxt_handle); 75827bd09bSSatish Balay } 76827bd09bSSatish Balay 777b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 7852f87cdaSBarry Smith PetscInt XXT_factor(xxt_ADT xxt_handle, /* prev. allocated xxt handle */ 7952f87cdaSBarry Smith PetscInt *local2global, /* global column mapping */ 8052f87cdaSBarry Smith PetscInt n, /* local num rows */ 8152f87cdaSBarry Smith PetscInt m, /* local num cols */ 82827bd09bSSatish Balay void *matvec, /* b_loc=A_local.x_loc */ 83827bd09bSSatish Balay void *grid_data /* grid data for matvec */ 84827bd09bSSatish Balay ) 85827bd09bSSatish Balay { 86b1c944f5SJed Brown PCTFS_comm_init(); 87827bd09bSSatish Balay check_handle(xxt_handle); 88827bd09bSSatish Balay 89827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 90b1c944f5SJed Brown if ((1<<(xxt_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); 91827bd09bSSatish Balay 92827bd09bSSatish Balay /* space for X info */ 93a501084fSBarry Smith xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info)); 94827bd09bSSatish Balay 95827bd09bSSatish Balay /* set up matvec handles */ 96827bd09bSSatish Balay xxt_handle->mvi = set_mvi(local2global, n, m, matvec, grid_data); 97827bd09bSSatish Balay 98827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 99827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 100827bd09bSSatish Balay xxt_handle->ns=0; 101827bd09bSSatish Balay 102827bd09bSSatish Balay /* determine separators and generate firing order - NB xxt info set here */ 103827bd09bSSatish Balay det_separators(xxt_handle); 104827bd09bSSatish Balay 105827bd09bSSatish Balay return(do_xxt_factor(xxt_handle)); 106827bd09bSSatish Balay } 107827bd09bSSatish Balay 1087b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 1098cda6cd7SBarry Smith PetscInt XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b) 110827bd09bSSatish Balay { 111827bd09bSSatish Balay 112b1c944f5SJed Brown PCTFS_comm_init(); 113827bd09bSSatish Balay check_handle(xxt_handle); 114827bd09bSSatish Balay 115827bd09bSSatish Balay /* need to copy b into x? */ 116827bd09bSSatish Balay if (b) 117*ca8e9878SJed Brown {PCTFS_rvec_copy(x,b,xxt_handle->mvi->n);} 118827bd09bSSatish Balay do_xxt_solve(xxt_handle,x); 119827bd09bSSatish Balay 120827bd09bSSatish Balay return(0); 121827bd09bSSatish Balay } 122827bd09bSSatish Balay 1237b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 12452f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle) 125827bd09bSSatish Balay { 126827bd09bSSatish Balay 127b1c944f5SJed Brown PCTFS_comm_init(); 128827bd09bSSatish Balay check_handle(xxt_handle); 129827bd09bSSatish Balay n_xxt_handles--; 130827bd09bSSatish Balay 131a501084fSBarry Smith free(xxt_handle->info->nsep); 132a501084fSBarry Smith free(xxt_handle->info->lnsep); 133a501084fSBarry Smith free(xxt_handle->info->fo); 134a501084fSBarry Smith free(xxt_handle->info->stages); 135a501084fSBarry Smith free(xxt_handle->info->solve_uu); 136a501084fSBarry Smith free(xxt_handle->info->solve_w); 137a501084fSBarry Smith free(xxt_handle->info->x); 138a501084fSBarry Smith free(xxt_handle->info->col_vals); 139a501084fSBarry Smith free(xxt_handle->info->col_sz); 140a501084fSBarry Smith free(xxt_handle->info->col_indices); 141a501084fSBarry Smith free(xxt_handle->info); 142a501084fSBarry Smith free(xxt_handle->mvi->local2global); 143*ca8e9878SJed Brown PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle); 144a501084fSBarry Smith free(xxt_handle->mvi); 145a501084fSBarry Smith free(xxt_handle); 146827bd09bSSatish Balay 147827bd09bSSatish Balay /* if the check fails we nuke */ 148a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 149827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 150827bd09bSSatish Balay return(0); 151827bd09bSSatish Balay } 152827bd09bSSatish Balay 1537b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 15452f87cdaSBarry Smith PetscInt XXT_stats(xxt_ADT xxt_handle) 155827bd09bSSatish Balay { 15652f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 15752f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 15852f87cdaSBarry Smith PetscInt vals[9], work[9]; 159a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 16071a0148aSBarry Smith PetscErrorCode ierr; 161827bd09bSSatish Balay 162b1c944f5SJed Brown PCTFS_comm_init(); 163827bd09bSSatish Balay check_handle(xxt_handle); 164827bd09bSSatish Balay 165827bd09bSSatish Balay /* if factorization not done there are no stats */ 166827bd09bSSatish Balay if (!xxt_handle->info||!xxt_handle->mvi) 167827bd09bSSatish Balay { 168b1c944f5SJed Brown if (!PCTFS_my_id) {ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");CHKERRQ(ierr);} 169827bd09bSSatish Balay return 1; 170827bd09bSSatish Balay } 171827bd09bSSatish Balay 172827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz; 173827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n; 174827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz; 175b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 176827bd09bSSatish Balay 177827bd09bSSatish Balay fvals[0]=fvals[1]=fvals[2] 178827bd09bSSatish Balay =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++; 179b1c944f5SJed Brown PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 180827bd09bSSatish Balay 181b1c944f5SJed Brown if (!PCTFS_my_id) 182827bd09bSSatish Balay { 183b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_nnz=%D\n",PCTFS_my_id,vals[0]);CHKERRQ(ierr); 184b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_nnz=%D\n",PCTFS_my_id,vals[1]);CHKERRQ(ierr); 185b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 186b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_nnz=%D\n",PCTFS_my_id,vals[2]);CHKERRQ(ierr); 187b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt C(2d) =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.5)));CHKERRQ(ierr); 188b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt C(3d) =%g\n",PCTFS_my_id,vals[2]/(pow(1.0*vals[5],1.6667)));CHKERRQ(ierr); 189b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_n =%D\n",PCTFS_my_id,vals[3]);CHKERRQ(ierr); 190b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_n =%D\n",PCTFS_my_id,vals[4]);CHKERRQ(ierr); 191b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);CHKERRQ(ierr); 192b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_n =%D\n",PCTFS_my_id,vals[5]);CHKERRQ(ierr); 193b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_buf=%D\n",PCTFS_my_id,vals[6]);CHKERRQ(ierr); 194b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_buf=%D\n",PCTFS_my_id,vals[7]);CHKERRQ(ierr); 195b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);CHKERRQ(ierr); 196b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_slv=%g\n",PCTFS_my_id,fvals[0]);CHKERRQ(ierr); 197b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_slv=%g\n",PCTFS_my_id,fvals[1]);CHKERRQ(ierr); 198b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 199827bd09bSSatish Balay } 200827bd09bSSatish Balay 201827bd09bSSatish Balay return(0); 202827bd09bSSatish Balay } 203827bd09bSSatish Balay 204827bd09bSSatish Balay /*************************************xxt.c************************************ 205827bd09bSSatish Balay 206827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 207827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 208827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 209827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 210827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 211827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 212*ca8e9878SJed Brown PCTFS_gs_init/gop). 213827bd09bSSatish Balay 214827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 215827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 216827bd09bSSatish Balay **************************************xxt.c***********************************/ 21752f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle) 218827bd09bSSatish Balay { 2197b1ae94cSBarry Smith return xxt_generate(xxt_handle); 220827bd09bSSatish Balay } 221827bd09bSSatish Balay 2227b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 22352f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle) 224827bd09bSSatish Balay { 22552f87cdaSBarry Smith PetscInt i,j,k,idex; 22652f87cdaSBarry Smith PetscInt dim, col; 227a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22852f87cdaSBarry Smith PetscInt *segs; 22952f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 23052f87cdaSBarry Smith PetscInt off, len; 231a501084fSBarry Smith PetscScalar *x_ptr; 23252f87cdaSBarry Smith PetscInt *iptr, flag; 23352f87cdaSBarry Smith PetscInt start=0, end, work; 23452f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 235*ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 23652f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 23752f87cdaSBarry Smith PetscInt a_n=xxt_handle->mvi->n; 23852f87cdaSBarry Smith PetscInt a_m=xxt_handle->mvi->m; 23952f87cdaSBarry Smith PetscInt *a_local2global=xxt_handle->mvi->local2global; 24052f87cdaSBarry Smith PetscInt level; 24152f87cdaSBarry Smith PetscInt xxt_nnz=0, xxt_max_nnz=0; 24252f87cdaSBarry Smith PetscInt n, m; 24352f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 244a501084fSBarry Smith PetscScalar **col_vals, *x; 24552f87cdaSBarry Smith PetscInt n_global; 24652f87cdaSBarry Smith PetscInt xxt_zero_nnz=0; 24752f87cdaSBarry Smith PetscInt xxt_zero_nnz_0=0; 24871a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 249a501084fSBarry Smith PetscScalar dm1 = -1.0; 250d1528f56SBarry Smith PetscErrorCode ierr; 251827bd09bSSatish Balay 252827bd09bSSatish Balay n=xxt_handle->mvi->n; 253827bd09bSSatish Balay nsep=xxt_handle->info->nsep; 254827bd09bSSatish Balay lnsep=xxt_handle->info->lnsep; 255827bd09bSSatish Balay fo=xxt_handle->info->fo; 256827bd09bSSatish Balay end=lnsep[0]; 257827bd09bSSatish Balay level=xxt_handle->level; 258*ca8e9878SJed Brown PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle; 259827bd09bSSatish Balay 260827bd09bSSatish Balay /* is there a null space? */ 261827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 262827bd09bSSatish Balay for (i=0, j=0; i<=level; i++) 263827bd09bSSatish Balay {j+=nsep[i];} 264827bd09bSSatish Balay 265827bd09bSSatish Balay m = j-xxt_handle->ns; 266827bd09bSSatish Balay if (m!=j) 26771a0148aSBarry Smith {ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);} 268827bd09bSSatish Balay 269827bd09bSSatish Balay /* get and initialize storage for x local */ 270827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 27152f87cdaSBarry Smith col_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 27252f87cdaSBarry Smith col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 273a501084fSBarry Smith col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *)); 274827bd09bSSatish Balay for (i=j=0; i<m; i++, j+=2) 275827bd09bSSatish Balay { 276827bd09bSSatish Balay col_indices[j]=col_indices[j+1]=col_sz[i]=-1; 277827bd09bSSatish Balay col_vals[i] = NULL; 278827bd09bSSatish Balay } 279827bd09bSSatish Balay col_indices[j]=-1; 280827bd09bSSatish Balay 281827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 282827bd09bSSatish Balay /* this looks like nsep[]=segments */ 28352f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 28452f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 285*ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 286*ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 287827bd09bSSatish Balay for (i=0; i<level; i++) 288827bd09bSSatish Balay {segs[i+1] += segs[i];} 289827bd09bSSatish Balay stages[0] = segs[0]; 290827bd09bSSatish Balay 291827bd09bSSatish Balay /* temporary vectors */ 292a501084fSBarry Smith u = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 293a501084fSBarry Smith z = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 294a501084fSBarry Smith v = (PetscScalar *) malloc(a_m*sizeof(PetscScalar)); 295a501084fSBarry Smith uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 296a501084fSBarry Smith w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 297827bd09bSSatish Balay 298827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 299827bd09bSSatish Balay for (i=1, j=0; i<=level; i++) 300827bd09bSSatish Balay {j+=nsep[i];} 301827bd09bSSatish Balay 302827bd09bSSatish Balay /* storage for sparse x values */ 303827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 304b1c944f5SJed Brown xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 305a501084fSBarry Smith x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar)); 306827bd09bSSatish Balay xxt_nnz = 0; 307827bd09bSSatish Balay 308827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 309827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 310827bd09bSSatish Balay for (dim=i=j=0;i<m;i++) 311827bd09bSSatish Balay { 312827bd09bSSatish Balay /* time to move to the next level? */ 313d4af0d30SBarry Smith while (i==segs[dim]) { 314e32f2f54SBarry Smith if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 315827bd09bSSatish Balay stages[dim++]=i; 316827bd09bSSatish Balay end+=lnsep[dim]; 317827bd09bSSatish Balay } 318827bd09bSSatish Balay stages[dim]=i; 319827bd09bSSatish Balay 320827bd09bSSatish Balay /* which column are we firing? */ 321827bd09bSSatish Balay /* i.e. set v_l */ 322827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 323827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 324b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 325827bd09bSSatish Balay 326827bd09bSSatish Balay /* shouldn't need this */ 327827bd09bSSatish Balay if (col==INT_MAX) 328827bd09bSSatish Balay { 329f1ed62a8SBarry Smith ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 330827bd09bSSatish Balay continue; 331827bd09bSSatish Balay } 332827bd09bSSatish Balay 333827bd09bSSatish Balay /* do I own it? I should */ 334*ca8e9878SJed Brown PCTFS_rvec_zero(v ,a_m); 335827bd09bSSatish Balay if (col==fo[start]) 336827bd09bSSatish Balay { 337827bd09bSSatish Balay start++; 338*ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_n); 339e7e72b3dSBarry Smith if (idex!=-1){ 340e7e72b3dSBarry Smith v[idex] = 1.0; j++; 341e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n"); 342827bd09bSSatish Balay } 343827bd09bSSatish Balay else 344827bd09bSSatish Balay { 345*ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_m); 346827bd09bSSatish Balay if (idex!=-1) 347827bd09bSSatish Balay {v[idex] = 1.0;} 348827bd09bSSatish Balay } 349827bd09bSSatish Balay 350827bd09bSSatish Balay /* perform u = A.v_l */ 351*ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 352827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 353827bd09bSSatish Balay 354827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 355827bd09bSSatish Balay /* technically only need to zero out first i entries */ 356827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 357*ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 358827bd09bSSatish Balay x_ptr=x; 359827bd09bSSatish Balay iptr = col_indices; 360827bd09bSSatish Balay for (k=0; k<i; k++) 361827bd09bSSatish Balay { 362827bd09bSSatish Balay off = *iptr++; 363827bd09bSSatish Balay len = *iptr++; 3640805154bSBarry Smith dlen = PetscBLASIntCast(len); 36571a0148aSBarry Smith uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1); 366827bd09bSSatish Balay x_ptr+=len; 367827bd09bSSatish Balay } 368827bd09bSSatish Balay 369827bd09bSSatish Balay 370827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 371b1c944f5SJed Brown PCTFS_ssgl_radd (uu, w, dim, stages); 372827bd09bSSatish Balay 373827bd09bSSatish Balay /* z = X.uu */ 374*ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 375827bd09bSSatish Balay x_ptr=x; 376827bd09bSSatish Balay iptr = col_indices; 377827bd09bSSatish Balay for (k=0; k<i; k++) 378827bd09bSSatish Balay { 379827bd09bSSatish Balay off = *iptr++; 380827bd09bSSatish Balay len = *iptr++; 3810805154bSBarry Smith dlen = PetscBLASIntCast(len); 38271a0148aSBarry Smith BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1); 383827bd09bSSatish Balay x_ptr+=len; 384827bd09bSSatish Balay } 385827bd09bSSatish Balay 386827bd09bSSatish Balay /* compute v_l = v_l - z */ 387*ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 3880805154bSBarry Smith dlen = PetscBLASIntCast(n); 38971a0148aSBarry Smith BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1); 390827bd09bSSatish Balay 391827bd09bSSatish Balay /* compute u_l = A.v_l */ 392827bd09bSSatish Balay if (a_n!=a_m) 393*ca8e9878SJed Brown {PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim);} 394*ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 395827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 396827bd09bSSatish Balay 397827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 3980805154bSBarry Smith dlen = PetscBLASIntCast(n); 39971a0148aSBarry Smith alpha = BLASdot_(&dlen,u,&i1,v,&i1); 400827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 401b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 402827bd09bSSatish Balay 4038f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 404827bd09bSSatish Balay 405827bd09bSSatish Balay /* check for small alpha */ 406827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 407c1235816SBarry Smith if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 408827bd09bSSatish Balay 409827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 410*ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 411827bd09bSSatish Balay 412827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 413827bd09bSSatish Balay flag = 1; 414827bd09bSSatish Balay off=len=0; 415827bd09bSSatish Balay for (k=0; k<n; k++) 416827bd09bSSatish Balay { 417827bd09bSSatish Balay if (v[k]!=0.0) 418827bd09bSSatish Balay { 419827bd09bSSatish Balay len=k; 420827bd09bSSatish Balay if (flag) 421827bd09bSSatish Balay {off=k; flag=0;} 422827bd09bSSatish Balay } 423827bd09bSSatish Balay } 424827bd09bSSatish Balay 425827bd09bSSatish Balay len -= (off-1); 426827bd09bSSatish Balay 427827bd09bSSatish Balay if (len>0) 428827bd09bSSatish Balay { 429827bd09bSSatish Balay if ((xxt_nnz+len)>xxt_max_nnz) 430827bd09bSSatish Balay { 431f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 432827bd09bSSatish Balay xxt_max_nnz *= 2; 433a501084fSBarry Smith x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar)); 434*ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xxt_nnz); 435a501084fSBarry Smith free(x); 436827bd09bSSatish Balay x = x_ptr; 437827bd09bSSatish Balay x_ptr+=xxt_nnz; 438827bd09bSSatish Balay } 439827bd09bSSatish Balay xxt_nnz += len; 440*ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 441827bd09bSSatish Balay 442827bd09bSSatish Balay /* keep track of number of zeros */ 443827bd09bSSatish Balay if (dim) 444827bd09bSSatish Balay { 445827bd09bSSatish Balay for (k=0; k<len; k++) 446827bd09bSSatish Balay { 447827bd09bSSatish Balay if (x_ptr[k]==0.0) 448827bd09bSSatish Balay {xxt_zero_nnz++;} 449827bd09bSSatish Balay } 450827bd09bSSatish Balay } 451827bd09bSSatish Balay else 452827bd09bSSatish Balay { 453827bd09bSSatish Balay for (k=0; k<len; k++) 454827bd09bSSatish Balay { 455827bd09bSSatish Balay if (x_ptr[k]==0.0) 456827bd09bSSatish Balay {xxt_zero_nnz_0++;} 457827bd09bSSatish Balay } 458827bd09bSSatish Balay } 459827bd09bSSatish Balay col_indices[2*i] = off; 460827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = len; 461827bd09bSSatish Balay col_vals[i] = x_ptr; 462827bd09bSSatish Balay } 463827bd09bSSatish Balay else 464827bd09bSSatish Balay { 465827bd09bSSatish Balay col_indices[2*i] = 0; 466827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = 0; 467827bd09bSSatish Balay col_vals[i] = x_ptr; 468827bd09bSSatish Balay } 469827bd09bSSatish Balay } 470827bd09bSSatish Balay 471827bd09bSSatish Balay /* close off stages for execution phase */ 472827bd09bSSatish Balay while (dim!=level) 473827bd09bSSatish Balay { 474827bd09bSSatish Balay stages[dim++]=i; 47571a0148aSBarry Smith ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 476827bd09bSSatish Balay } 477827bd09bSSatish Balay stages[dim]=i; 478827bd09bSSatish Balay 479827bd09bSSatish Balay xxt_handle->info->n=xxt_handle->mvi->n; 480827bd09bSSatish Balay xxt_handle->info->m=m; 481827bd09bSSatish Balay xxt_handle->info->nnz=xxt_nnz; 482827bd09bSSatish Balay xxt_handle->info->max_nnz=xxt_max_nnz; 483827bd09bSSatish Balay xxt_handle->info->msg_buf_sz=stages[level]-stages[0]; 484a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 485a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 486827bd09bSSatish Balay xxt_handle->info->x=x; 487827bd09bSSatish Balay xxt_handle->info->col_vals=col_vals; 488827bd09bSSatish Balay xxt_handle->info->col_sz=col_sz; 489827bd09bSSatish Balay xxt_handle->info->col_indices=col_indices; 490827bd09bSSatish Balay xxt_handle->info->stages=stages; 491827bd09bSSatish Balay xxt_handle->info->nsolves=0; 492827bd09bSSatish Balay xxt_handle->info->tot_solve_time=0.0; 493827bd09bSSatish Balay 494a501084fSBarry Smith free(segs); 495a501084fSBarry Smith free(u); 496a501084fSBarry Smith free(v); 497a501084fSBarry Smith free(uu); 498a501084fSBarry Smith free(z); 499a501084fSBarry Smith free(w); 500827bd09bSSatish Balay 501827bd09bSSatish Balay return(0); 502827bd09bSSatish Balay } 503827bd09bSSatish Balay 5047b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5050924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) 506827bd09bSSatish Balay { 50752f87cdaSBarry Smith PetscInt off, len, *iptr; 50852f87cdaSBarry Smith PetscInt level =xxt_handle->level; 50952f87cdaSBarry Smith PetscInt n =xxt_handle->info->n; 51052f87cdaSBarry Smith PetscInt m =xxt_handle->info->m; 51152f87cdaSBarry Smith PetscInt *stages =xxt_handle->info->stages; 51252f87cdaSBarry Smith PetscInt *col_indices=xxt_handle->info->col_indices; 513a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 514a501084fSBarry Smith PetscScalar *solve_uu=xxt_handle->info->solve_uu; 515a501084fSBarry Smith PetscScalar *solve_w =xxt_handle->info->solve_w; 516a501084fSBarry Smith PetscScalar *x =xxt_handle->info->x; 51771a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 518827bd09bSSatish Balay 5190924e98cSBarry Smith PetscFunctionBegin; 520827bd09bSSatish Balay uu_ptr=solve_uu; 521*ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 522827bd09bSSatish Balay 523827bd09bSSatish Balay /* x = X.Y^T.b */ 524827bd09bSSatish Balay /* uu = Y^T.b */ 525827bd09bSSatish Balay for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) 526827bd09bSSatish Balay { 527827bd09bSSatish Balay off=*iptr++; len=*iptr++; 5280805154bSBarry Smith dlen = PetscBLASIntCast(len); 52971a0148aSBarry Smith *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1); 530827bd09bSSatish Balay } 531827bd09bSSatish Balay 532827bd09bSSatish Balay /* comunication of beta */ 533827bd09bSSatish Balay uu_ptr=solve_uu; 534b1c944f5SJed Brown if (level) {PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);} 535827bd09bSSatish Balay 536*ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 537827bd09bSSatish Balay 538827bd09bSSatish Balay /* x = X.uu */ 539827bd09bSSatish Balay for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) 540827bd09bSSatish Balay { 541827bd09bSSatish Balay off=*iptr++; len=*iptr++; 5420805154bSBarry Smith dlen = PetscBLASIntCast(len); 54371a0148aSBarry Smith BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1); 544827bd09bSSatish Balay } 5450924e98cSBarry Smith PetscFunctionReturn(0); 546827bd09bSSatish Balay } 547827bd09bSSatish Balay 5487b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5490924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle) 550827bd09bSSatish Balay { 55152f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 552827bd09bSSatish Balay 5530924e98cSBarry Smith PetscFunctionBegin; 554e32f2f54SBarry Smith if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle); 555827bd09bSSatish Balay 556827bd09bSSatish Balay vals[0]=vals[1]=xxt_handle->id; 557b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 558e32f2f54SBarry Smith if ((vals[0]!=vals[1])||(xxt_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], xxt_handle->id); 5590924e98cSBarry Smith PetscFunctionReturn(0); 560827bd09bSSatish Balay } 561827bd09bSSatish Balay 5627b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5630924e98cSBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle) 564827bd09bSSatish Balay { 56552f87cdaSBarry Smith PetscInt i, ct, id; 56652f87cdaSBarry Smith PetscInt mask, edge, *iptr; 56752f87cdaSBarry Smith PetscInt *dir, *used; 56852f87cdaSBarry Smith PetscInt sum[4], w[4]; 569a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 57052f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 571a501084fSBarry Smith PetscScalar *lhs, *rhs; 57252f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 573*ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle; 57452f87cdaSBarry Smith PetscInt *local2global=xxt_handle->mvi->local2global; 57552f87cdaSBarry Smith PetscInt n=xxt_handle->mvi->n; 57652f87cdaSBarry Smith PetscInt m=xxt_handle->mvi->m; 57752f87cdaSBarry Smith PetscInt level=xxt_handle->level; 578ab824b78SBarry Smith PetscInt shared=0; 579827bd09bSSatish Balay 5800924e98cSBarry Smith PetscFunctionBegin; 58152f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 58252f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 58352f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 58452f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 58552f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 586827bd09bSSatish Balay 587*ca8e9878SJed Brown PCTFS_ivec_zero(dir ,level+1); 588*ca8e9878SJed Brown PCTFS_ivec_zero(nsep ,level+1); 589*ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 590*ca8e9878SJed Brown PCTFS_ivec_set (fo ,-1,n+1); 591*ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 592827bd09bSSatish Balay 5938cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 5948cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 595827bd09bSSatish Balay 596827bd09bSSatish Balay /* determine the # of unique dof */ 597*ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 598*ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 599*ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 600*ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 601827bd09bSSatish Balay for (ct=i=0;i<n;i++) 602827bd09bSSatish Balay { 603827bd09bSSatish Balay if (lhs[i]!=0.0) 604827bd09bSSatish Balay {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];} 605827bd09bSSatish Balay } 606b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 607827bd09bSSatish Balay rsum[0]+=0.1; 608827bd09bSSatish Balay rsum[1]+=0.1; 609827bd09bSSatish Balay 610827bd09bSSatish Balay if (fabs(rsum[0]-rsum[1])>EPS) 611ab824b78SBarry Smith {shared=1;} 612827bd09bSSatish Balay 61352f87cdaSBarry Smith xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0]; 61452f87cdaSBarry Smith xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0]; 615827bd09bSSatish Balay 616827bd09bSSatish Balay /* determine separator sets top down */ 617827bd09bSSatish Balay if (shared) 618827bd09bSSatish Balay { 619b1c944f5SJed Brown for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 620827bd09bSSatish Balay { 621827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 622*ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 623*ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 624827bd09bSSatish Balay 625827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 626*ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 627*ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 628827bd09bSSatish Balay 629827bd09bSSatish Balay for (i=0;i<n;i++) 630827bd09bSSatish Balay { 631827bd09bSSatish Balay if (id< mask) 632827bd09bSSatish Balay { 633827bd09bSSatish Balay if (lhs[i]!=0.0) 634827bd09bSSatish Balay {lhs[i]=1.0;} 635827bd09bSSatish Balay } 636827bd09bSSatish Balay if (id>=mask) 637827bd09bSSatish Balay { 638827bd09bSSatish Balay if (rhs[i]!=0.0) 639827bd09bSSatish Balay {rhs[i]=1.0;} 640827bd09bSSatish Balay } 641827bd09bSSatish Balay } 642827bd09bSSatish Balay 643827bd09bSSatish Balay if (id< mask) 644*ca8e9878SJed Brown {PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1);} 645827bd09bSSatish Balay else 646*ca8e9878SJed Brown {PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1);} 647827bd09bSSatish Balay 648827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 649*ca8e9878SJed Brown PCTFS_rvec_zero(rsum,4); 650*ca8e9878SJed Brown for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) 651827bd09bSSatish Balay { 652827bd09bSSatish Balay if (!used[i]) 653827bd09bSSatish Balay { 654827bd09bSSatish Balay /* number of unmarked dofs on node */ 655827bd09bSSatish Balay ct++; 656827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 657827bd09bSSatish Balay if (id< mask) 658827bd09bSSatish Balay { 659827bd09bSSatish Balay if (lhs[i]!=0.0) 660827bd09bSSatish Balay {sum[0]++; rsum[0]+=1.0/lhs[i];} 661827bd09bSSatish Balay } 662827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 663827bd09bSSatish Balay if (id>=mask) 664827bd09bSSatish Balay { 665827bd09bSSatish Balay if (rhs[i]!=0.0) 666827bd09bSSatish Balay {sum[1]++; rsum[1]+=1.0/rhs[i];} 667827bd09bSSatish Balay } 668827bd09bSSatish Balay } 669827bd09bSSatish Balay } 670827bd09bSSatish Balay 671827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 672827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 673827bd09bSSatish Balay (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct); 674b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 675b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,4,op,edge); 676827bd09bSSatish Balay rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1; 677827bd09bSSatish Balay 678827bd09bSSatish Balay if (id<mask) 679827bd09bSSatish Balay { 680827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 681827bd09bSSatish Balay for (ct=i=0;i<n;i++) 682827bd09bSSatish Balay { 683827bd09bSSatish Balay if ((!used[i])&&(lhs[i]!=0.0)) 684827bd09bSSatish Balay { 685827bd09bSSatish Balay ct++; nfo++; 686827bd09bSSatish Balay 687c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 688827bd09bSSatish Balay 689827bd09bSSatish Balay *--iptr = local2global[i]; 690827bd09bSSatish Balay used[i]=edge; 691827bd09bSSatish Balay } 692827bd09bSSatish Balay } 693*ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 694827bd09bSSatish Balay 695827bd09bSSatish Balay lnsep[edge]=ct; 69652f87cdaSBarry Smith nsep[edge]=(PetscInt) rsum[0]; 697827bd09bSSatish Balay dir [edge]=LEFT; 698827bd09bSSatish Balay } 699827bd09bSSatish Balay 700827bd09bSSatish Balay if (id>=mask) 701827bd09bSSatish Balay { 702827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 703827bd09bSSatish Balay for (ct=i=0;i<n;i++) 704827bd09bSSatish Balay { 705827bd09bSSatish Balay if ((!used[i])&&(rhs[i]!=0.0)) 706827bd09bSSatish Balay { 707827bd09bSSatish Balay ct++; nfo++; 708827bd09bSSatish Balay 709c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 710827bd09bSSatish Balay 711827bd09bSSatish Balay *--iptr = local2global[i]; 712827bd09bSSatish Balay used[i]=edge; 713827bd09bSSatish Balay } 714827bd09bSSatish Balay } 715*ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 716827bd09bSSatish Balay 717827bd09bSSatish Balay lnsep[edge]=ct; 71852f87cdaSBarry Smith nsep[edge]= (PetscInt) rsum[1]; 719827bd09bSSatish Balay dir [edge]=RIGHT; 720827bd09bSSatish Balay } 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 */ 726827bd09bSSatish Balay if (id>=mask) 727827bd09bSSatish Balay {id-=mask;} 728827bd09bSSatish Balay } 729827bd09bSSatish Balay } 730827bd09bSSatish Balay else 731827bd09bSSatish Balay { 732b1c944f5SJed Brown for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) 733827bd09bSSatish Balay { 734827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 735*ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 736*ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 737827bd09bSSatish Balay 738827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 739*ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 740*ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 741827bd09bSSatish Balay 742827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 743*ca8e9878SJed Brown for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) 744827bd09bSSatish Balay { 745827bd09bSSatish Balay if (!used[i]) 746827bd09bSSatish Balay { 747827bd09bSSatish Balay /* number of unmarked dofs on node */ 748827bd09bSSatish Balay ct++; 749827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 750827bd09bSSatish Balay if ((id< mask)&&(lhs[i]!=0.0)) {sum[0]++;} 751827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 752827bd09bSSatish Balay if ((id>=mask)&&(rhs[i]!=0.0)) {sum[1]++;} 753827bd09bSSatish Balay } 754827bd09bSSatish Balay } 755827bd09bSSatish Balay 756827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 757827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 758b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 759827bd09bSSatish Balay 760827bd09bSSatish Balay /* lhs hc wins */ 761827bd09bSSatish Balay if (sum[2]>=sum[3]) 762827bd09bSSatish Balay { 763827bd09bSSatish Balay if (id<mask) 764827bd09bSSatish Balay { 765827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 766827bd09bSSatish Balay for (ct=i=0;i<n;i++) 767827bd09bSSatish Balay { 768827bd09bSSatish Balay if ((!used[i])&&(lhs[i]!=0.0)) 769827bd09bSSatish Balay { 770827bd09bSSatish Balay ct++; nfo++; 771827bd09bSSatish Balay *--iptr = local2global[i]; 772827bd09bSSatish Balay used[i]=edge; 773827bd09bSSatish Balay } 774827bd09bSSatish Balay } 775*ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 776827bd09bSSatish Balay lnsep[edge]=ct; 777827bd09bSSatish Balay } 778827bd09bSSatish Balay nsep[edge]=sum[0]; 779827bd09bSSatish Balay dir [edge]=LEFT; 780827bd09bSSatish Balay } 781827bd09bSSatish Balay /* rhs hc wins */ 782827bd09bSSatish Balay else 783827bd09bSSatish Balay { 784827bd09bSSatish Balay if (id>=mask) 785827bd09bSSatish Balay { 786827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 787827bd09bSSatish Balay for (ct=i=0;i<n;i++) 788827bd09bSSatish Balay { 789827bd09bSSatish Balay if ((!used[i])&&(rhs[i]!=0.0)) 790827bd09bSSatish Balay { 791827bd09bSSatish Balay ct++; nfo++; 792827bd09bSSatish Balay *--iptr = local2global[i]; 793827bd09bSSatish Balay used[i]=edge; 794827bd09bSSatish Balay } 795827bd09bSSatish Balay } 796*ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 797827bd09bSSatish Balay lnsep[edge]=ct; 798827bd09bSSatish Balay } 799827bd09bSSatish Balay nsep[edge]=sum[1]; 800827bd09bSSatish Balay dir [edge]=RIGHT; 801827bd09bSSatish Balay } 802827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 803827bd09bSSatish Balay /* do we need full set of separators for this? */ 804827bd09bSSatish Balay 805827bd09bSSatish Balay /* fold rhs hc into lower */ 806827bd09bSSatish Balay if (id>=mask) 807827bd09bSSatish Balay {id-=mask;} 808827bd09bSSatish Balay } 809827bd09bSSatish Balay } 810827bd09bSSatish Balay 811827bd09bSSatish Balay 812827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 813827bd09bSSatish Balay for (ct=i=0;i<n;i++) 814827bd09bSSatish Balay { 815827bd09bSSatish Balay if (!used[i]) 816827bd09bSSatish Balay { 817827bd09bSSatish Balay ct++; nfo++; 818827bd09bSSatish Balay *--iptr = local2global[i]; 819827bd09bSSatish Balay used[i]=edge; 820827bd09bSSatish Balay } 821827bd09bSSatish Balay } 822*ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 823827bd09bSSatish Balay lnsep[edge]=ct; 824827bd09bSSatish Balay nsep [edge]=ct; 825827bd09bSSatish Balay dir [edge]=LEFT; 826827bd09bSSatish Balay 827827bd09bSSatish Balay xxt_handle->info->nsep=nsep; 828827bd09bSSatish Balay xxt_handle->info->lnsep=lnsep; 829827bd09bSSatish Balay xxt_handle->info->fo=fo; 830827bd09bSSatish Balay xxt_handle->info->nfo=nfo; 831827bd09bSSatish Balay 832a501084fSBarry Smith free(dir); 833a501084fSBarry Smith free(lhs); 834a501084fSBarry Smith free(rhs); 835a501084fSBarry Smith free(used); 8360924e98cSBarry Smith PetscFunctionReturn(0); 837827bd09bSSatish Balay } 838827bd09bSSatish Balay 8397b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 84052f87cdaSBarry Smith static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, void *matvec, void *grid_data) 841827bd09bSSatish Balay { 842827bd09bSSatish Balay mv_info *mvi; 843827bd09bSSatish Balay 844827bd09bSSatish Balay 845a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 846827bd09bSSatish Balay mvi->n=n; 847827bd09bSSatish Balay mvi->m=m; 848827bd09bSSatish Balay mvi->n_global=-1; 849827bd09bSSatish Balay mvi->m_global=-1; 85052f87cdaSBarry Smith mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt)); 851*ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 852827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 853a501084fSBarry Smith mvi->matvec=(PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec; 854827bd09bSSatish Balay mvi->grid_data=grid_data; 855827bd09bSSatish Balay 856827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 857*ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 858827bd09bSSatish Balay 859827bd09bSSatish Balay return(mvi); 860827bd09bSSatish Balay } 861827bd09bSSatish Balay 8627b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 8630924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 864827bd09bSSatish Balay { 8650924e98cSBarry Smith PetscFunctionBegin; 866827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 8670924e98cSBarry Smith PetscFunctionReturn(0); 868827bd09bSSatish Balay } 869827bd09bSSatish Balay 870827bd09bSSatish Balay 871827bd09bSSatish Balay 872