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; 38ca8e9878SJed 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); 59*38dcbb09SBarry Smith static PetscErrorCode xxt_generate(xxt_ADT xxt_handle); 60*38dcbb09SBarry Smith static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle); 615c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global, PetscInt n, PetscInt m, PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*), 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***********************************/ 78*38dcbb09SBarry Smith PetscErrorCode 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 */ 825c8f6a95SKarl Rupp PetscErrorCode (*matvec)(void*,PetscScalar*,PetscScalar*), /* b_loc=A_local.x_loc */ 831147fc2aSKarl Rupp void *grid_data) /* grid data for matvec */ 84827bd09bSSatish Balay { 85b1c944f5SJed Brown PCTFS_comm_init(); 86827bd09bSSatish Balay check_handle(xxt_handle); 87827bd09bSSatish Balay 88827bd09bSSatish Balay /* only 2^k for now and all nodes participating */ 89b1c944f5SJed 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); 90827bd09bSSatish Balay 91827bd09bSSatish Balay /* space for X info */ 92a501084fSBarry Smith xxt_handle->info = (xxt_info*)malloc(sizeof(xxt_info)); 93827bd09bSSatish Balay 94827bd09bSSatish Balay /* set up matvec handles */ 955c8f6a95SKarl Rupp xxt_handle->mvi = set_mvi(local2global, n, m, (PetscErrorCode (*)(mv_info*,PetscScalar*,PetscScalar*))matvec, grid_data); 96827bd09bSSatish Balay 97827bd09bSSatish Balay /* matrix is assumed to be of full rank */ 98827bd09bSSatish Balay /* LATER we can reset to indicate rank def. */ 99827bd09bSSatish Balay xxt_handle->ns=0; 100827bd09bSSatish Balay 101827bd09bSSatish Balay /* determine separators and generate firing order - NB xxt info set here */ 102827bd09bSSatish Balay det_separators(xxt_handle); 103827bd09bSSatish Balay 104827bd09bSSatish Balay return(do_xxt_factor(xxt_handle)); 105827bd09bSSatish Balay } 106827bd09bSSatish Balay 1077b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 108*38dcbb09SBarry Smith PetscErrorCode XXT_solve(xxt_ADT xxt_handle, PetscScalar *x, PetscScalar *b) 109827bd09bSSatish Balay { 110827bd09bSSatish Balay 111b1c944f5SJed Brown PCTFS_comm_init(); 112827bd09bSSatish Balay check_handle(xxt_handle); 113827bd09bSSatish Balay 114827bd09bSSatish Balay /* need to copy b into x? */ 1152fa5cd67SKarl Rupp if (b) PCTFS_rvec_copy(x,b,xxt_handle->mvi->n); 116*38dcbb09SBarry Smith return do_xxt_solve(xxt_handle,x); 117827bd09bSSatish Balay } 118827bd09bSSatish Balay 1197b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 12052f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle) 121827bd09bSSatish Balay { 122827bd09bSSatish Balay 123b1c944f5SJed Brown PCTFS_comm_init(); 124827bd09bSSatish Balay check_handle(xxt_handle); 125827bd09bSSatish Balay n_xxt_handles--; 126827bd09bSSatish Balay 127a501084fSBarry Smith free(xxt_handle->info->nsep); 128a501084fSBarry Smith free(xxt_handle->info->lnsep); 129a501084fSBarry Smith free(xxt_handle->info->fo); 130a501084fSBarry Smith free(xxt_handle->info->stages); 131a501084fSBarry Smith free(xxt_handle->info->solve_uu); 132a501084fSBarry Smith free(xxt_handle->info->solve_w); 133a501084fSBarry Smith free(xxt_handle->info->x); 134a501084fSBarry Smith free(xxt_handle->info->col_vals); 135a501084fSBarry Smith free(xxt_handle->info->col_sz); 136a501084fSBarry Smith free(xxt_handle->info->col_indices); 137a501084fSBarry Smith free(xxt_handle->info); 138a501084fSBarry Smith free(xxt_handle->mvi->local2global); 139ca8e9878SJed Brown PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle); 140a501084fSBarry Smith free(xxt_handle->mvi); 141a501084fSBarry Smith free(xxt_handle); 142827bd09bSSatish Balay 143827bd09bSSatish Balay /* if the check fails we nuke */ 144a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 145827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 146827bd09bSSatish Balay return(0); 147827bd09bSSatish Balay } 148827bd09bSSatish Balay 1497b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 150*38dcbb09SBarry Smith PetscErrorCode XXT_stats(xxt_ADT xxt_handle) 151827bd09bSSatish Balay { 15252f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 15352f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 15452f87cdaSBarry Smith PetscInt vals[9], work[9]; 155a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 15671a0148aSBarry Smith PetscErrorCode ierr; 157827bd09bSSatish Balay 158b1c944f5SJed Brown PCTFS_comm_init(); 159827bd09bSSatish Balay check_handle(xxt_handle); 160827bd09bSSatish Balay 161827bd09bSSatish Balay /* if factorization not done there are no stats */ 162db4deed7SKarl Rupp if (!xxt_handle->info||!xxt_handle->mvi) { 163b1c944f5SJed Brown if (!PCTFS_my_id) { ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");CHKERRQ(ierr); } 164827bd09bSSatish Balay return 1; 165827bd09bSSatish Balay } 166827bd09bSSatish Balay 167827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz; 168827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n; 169827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz; 170b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 171827bd09bSSatish Balay 1722fa5cd67SKarl Rupp fvals[0]=fvals[1]=fvals[2] =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++; 173b1c944f5SJed Brown PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 174827bd09bSSatish Balay 175db4deed7SKarl Rupp if (!PCTFS_my_id) { 176b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_nnz=%D\n",PCTFS_my_id,vals[0]);CHKERRQ(ierr); 177b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_nnz=%D\n",PCTFS_my_id,vals[1]);CHKERRQ(ierr); 178b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 179b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_nnz=%D\n",PCTFS_my_id,vals[2]);CHKERRQ(ierr); 18077b4d14cSPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt C(2d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.5)));CHKERRQ(ierr); 18177b4d14cSPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: xxt C(3d) =%g\n",PCTFS_my_id,vals[2]/(PetscPowReal(1.0*vals[5],1.6667)));CHKERRQ(ierr); 182b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_n =%D\n",PCTFS_my_id,vals[3]);CHKERRQ(ierr); 183b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_n =%D\n",PCTFS_my_id,vals[4]);CHKERRQ(ierr); 184b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);CHKERRQ(ierr); 185b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_n =%D\n",PCTFS_my_id,vals[5]);CHKERRQ(ierr); 186b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_buf=%D\n",PCTFS_my_id,vals[6]);CHKERRQ(ierr); 187b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_buf=%D\n",PCTFS_my_id,vals[7]);CHKERRQ(ierr); 188b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);CHKERRQ(ierr); 189b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_slv=%g\n",PCTFS_my_id,fvals[0]);CHKERRQ(ierr); 190b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_slv=%g\n",PCTFS_my_id,fvals[1]);CHKERRQ(ierr); 191b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 192827bd09bSSatish Balay } 193827bd09bSSatish Balay 194827bd09bSSatish Balay return(0); 195827bd09bSSatish Balay } 196827bd09bSSatish Balay 197827bd09bSSatish Balay /*************************************xxt.c************************************ 198827bd09bSSatish Balay 199827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 200827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 201827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 202827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 203827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 204827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 205ca8e9878SJed Brown PCTFS_gs_init/gop). 206827bd09bSSatish Balay 207827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 208827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 209827bd09bSSatish Balay **************************************xxt.c***********************************/ 210*38dcbb09SBarry Smith static PetscErrorCode do_xxt_factor(xxt_ADT xxt_handle) 211827bd09bSSatish Balay { 2127b1ae94cSBarry Smith return xxt_generate(xxt_handle); 213827bd09bSSatish Balay } 214827bd09bSSatish Balay 2157b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 216*38dcbb09SBarry Smith static PetscErrorCode xxt_generate(xxt_ADT xxt_handle) 217827bd09bSSatish Balay { 21852f87cdaSBarry Smith PetscInt i,j,k,idex; 21952f87cdaSBarry Smith PetscInt dim, col; 220a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22152f87cdaSBarry Smith PetscInt *segs; 22252f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 22352f87cdaSBarry Smith PetscInt off, len; 224a501084fSBarry Smith PetscScalar *x_ptr; 22552f87cdaSBarry Smith PetscInt *iptr, flag; 22652f87cdaSBarry Smith PetscInt start =0, end, work; 22752f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 228ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 22952f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 23052f87cdaSBarry Smith PetscInt a_n =xxt_handle->mvi->n; 23152f87cdaSBarry Smith PetscInt a_m =xxt_handle->mvi->m; 23252f87cdaSBarry Smith PetscInt *a_local2global=xxt_handle->mvi->local2global; 23352f87cdaSBarry Smith PetscInt level; 23452f87cdaSBarry Smith PetscInt xxt_nnz=0, xxt_max_nnz=0; 23552f87cdaSBarry Smith PetscInt n, m; 23652f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 237a501084fSBarry Smith PetscScalar **col_vals, *x; 23852f87cdaSBarry Smith PetscInt n_global; 23952f87cdaSBarry Smith PetscInt xxt_zero_nnz =0; 24052f87cdaSBarry Smith PetscInt xxt_zero_nnz_0=0; 24171a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 242a501084fSBarry Smith PetscScalar dm1 = -1.0; 243d1528f56SBarry Smith PetscErrorCode ierr; 244827bd09bSSatish Balay 245827bd09bSSatish Balay n = xxt_handle->mvi->n; 246827bd09bSSatish Balay nsep = xxt_handle->info->nsep; 247827bd09bSSatish Balay lnsep = xxt_handle->info->lnsep; 248827bd09bSSatish Balay fo = xxt_handle->info->fo; 249827bd09bSSatish Balay end = lnsep[0]; 250827bd09bSSatish Balay level = xxt_handle->level; 251ca8e9878SJed Brown PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 252827bd09bSSatish Balay 253827bd09bSSatish Balay /* is there a null space? */ 254827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 2552fa5cd67SKarl Rupp for (i=0, j=0; i<=level; i++) j+=nsep[i]; 256827bd09bSSatish Balay 257827bd09bSSatish Balay m = j-xxt_handle->ns; 25822d28d08SBarry Smith if (m!=j) { 25922d28d08SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);CHKERRQ(ierr); 26022d28d08SBarry Smith } 261827bd09bSSatish Balay 262827bd09bSSatish Balay /* get and initialize storage for x local */ 263827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 26452f87cdaSBarry Smith col_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 26552f87cdaSBarry Smith col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 266a501084fSBarry Smith col_vals = (PetscScalar**) malloc(m*sizeof(PetscScalar*)); 267db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 268827bd09bSSatish Balay col_indices[j]=col_indices[j+1]=col_sz[i]=-1; 269827bd09bSSatish Balay col_vals[i] = NULL; 270827bd09bSSatish Balay } 271827bd09bSSatish Balay col_indices[j]=-1; 272827bd09bSSatish Balay 273827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 274827bd09bSSatish Balay /* this looks like nsep[]=segments */ 27552f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 27652f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 277ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 278ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 2792fa5cd67SKarl Rupp for (i=0; i<level; i++) segs[i+1] += segs[i]; 280827bd09bSSatish Balay stages[0] = segs[0]; 281827bd09bSSatish Balay 282827bd09bSSatish Balay /* temporary vectors */ 283a501084fSBarry Smith u = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 284a501084fSBarry Smith z = (PetscScalar*) malloc(n*sizeof(PetscScalar)); 285a501084fSBarry Smith v = (PetscScalar*) malloc(a_m*sizeof(PetscScalar)); 286a501084fSBarry Smith uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 287a501084fSBarry Smith w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 288827bd09bSSatish Balay 289827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 2902fa5cd67SKarl Rupp for (i=1, j=0; i<=level; i++) j+=nsep[i]; 291827bd09bSSatish Balay 292827bd09bSSatish Balay /* storage for sparse x values */ 293827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 29485ec1a3cSBarry Smith xxt_max_nnz = (PetscInt)(2.5*PetscPowReal(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 295a501084fSBarry Smith x = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar)); 296827bd09bSSatish Balay xxt_nnz = 0; 297827bd09bSSatish Balay 298827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 299827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 3002fa5cd67SKarl Rupp for (dim=i=j=0; i<m; i++) { 301827bd09bSSatish Balay /* time to move to the next level? */ 302d4af0d30SBarry Smith while (i==segs[dim]) { 303e32f2f54SBarry Smith if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 304827bd09bSSatish Balay stages[dim++]=i; 305827bd09bSSatish Balay end +=lnsep[dim]; 306827bd09bSSatish Balay } 307827bd09bSSatish Balay stages[dim]=i; 308827bd09bSSatish Balay 309827bd09bSSatish Balay /* which column are we firing? */ 310827bd09bSSatish Balay /* i.e. set v_l */ 311827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 312827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 313b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 314827bd09bSSatish Balay 315827bd09bSSatish Balay /* shouldn't need this */ 316db4deed7SKarl Rupp if (col==INT_MAX) { 317f1ed62a8SBarry Smith ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 318827bd09bSSatish Balay continue; 319827bd09bSSatish Balay } 320827bd09bSSatish Balay 321827bd09bSSatish Balay /* do I own it? I should */ 322ca8e9878SJed Brown PCTFS_rvec_zero(v,a_m); 323db4deed7SKarl Rupp if (col==fo[start]) { 324827bd09bSSatish Balay start++; 325ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_n); 326e7e72b3dSBarry Smith if (idex!=-1) { 327e7e72b3dSBarry Smith v[idex] = 1.0; j++; 328e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n"); 329db4deed7SKarl Rupp } else { 330ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_m); 3312fa5cd67SKarl Rupp if (idex!=-1) v[idex] = 1.0; 332827bd09bSSatish Balay } 333827bd09bSSatish Balay 334827bd09bSSatish Balay /* perform u = A.v_l */ 335ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 336827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 337827bd09bSSatish Balay 338827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 339827bd09bSSatish Balay /* technically only need to zero out first i entries */ 340827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 341ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 342827bd09bSSatish Balay x_ptr=x; 343827bd09bSSatish Balay iptr = col_indices; 344db4deed7SKarl Rupp for (k=0; k<i; k++) { 345827bd09bSSatish Balay off = *iptr++; 346827bd09bSSatish Balay len = *iptr++; 347c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3488b83055fSJed Brown PetscStackCallBLAS("BLASdot",uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1)); 349827bd09bSSatish Balay x_ptr+=len; 350827bd09bSSatish Balay } 351827bd09bSSatish Balay 352827bd09bSSatish Balay 353827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 354a00740a5SBarry Smith ierr = PCTFS_ssgl_radd (uu, w, dim, stages);CHKERRQ(ierr); 355827bd09bSSatish Balay 356827bd09bSSatish Balay /* z = X.uu */ 357ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 358827bd09bSSatish Balay x_ptr=x; 359827bd09bSSatish Balay iptr = col_indices; 360db4deed7SKarl Rupp for (k=0; k<i; k++) { 361827bd09bSSatish Balay off = *iptr++; 362827bd09bSSatish Balay len = *iptr++; 363c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 3648b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1)); 365827bd09bSSatish Balay x_ptr+=len; 366827bd09bSSatish Balay } 367827bd09bSSatish Balay 368827bd09bSSatish Balay /* compute v_l = v_l - z */ 369ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 370c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 3718b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1)); 372827bd09bSSatish Balay 373827bd09bSSatish Balay /* compute u_l = A.v_l */ 3742fa5cd67SKarl Rupp if (a_n!=a_m) PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); 375ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 376827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 377827bd09bSSatish Balay 378827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 379c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&dlen);CHKERRQ(ierr); 3808b83055fSJed Brown PetscStackCallBLAS("BLASdot",alpha = BLASdot_(&dlen,u,&i1,v,&i1)); 381827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 382b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 383827bd09bSSatish Balay 3848f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 385827bd09bSSatish Balay 386827bd09bSSatish Balay /* check for small alpha */ 387827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 38877b4d14cSPeter Brune if (PetscAbsScalar(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 389827bd09bSSatish Balay 390827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 391ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 392827bd09bSSatish Balay 393827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 394827bd09bSSatish Balay flag = 1; 395827bd09bSSatish Balay off=len=0; 396db4deed7SKarl Rupp for (k=0; k<n; k++) { 397db4deed7SKarl Rupp if (v[k]!=0.0) { 398827bd09bSSatish Balay len=k; 399db4deed7SKarl Rupp if (flag) { off=k; flag=0; } 400827bd09bSSatish Balay } 401827bd09bSSatish Balay } 402827bd09bSSatish Balay 403827bd09bSSatish Balay len -= (off-1); 404827bd09bSSatish Balay 405db4deed7SKarl Rupp if (len>0) { 406db4deed7SKarl Rupp if ((xxt_nnz+len)>xxt_max_nnz) { 407f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 408827bd09bSSatish Balay xxt_max_nnz *= 2; 409a501084fSBarry Smith x_ptr = (PetscScalar*) malloc(xxt_max_nnz*sizeof(PetscScalar)); 410ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xxt_nnz); 411a501084fSBarry Smith free(x); 412827bd09bSSatish Balay x = x_ptr; 413827bd09bSSatish Balay x_ptr+=xxt_nnz; 414827bd09bSSatish Balay } 415827bd09bSSatish Balay xxt_nnz += len; 416ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 417827bd09bSSatish Balay 418827bd09bSSatish Balay /* keep track of number of zeros */ 419db4deed7SKarl Rupp if (dim) { 420db4deed7SKarl Rupp for (k=0; k<len; k++) { 4212fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xxt_zero_nnz++; 422827bd09bSSatish Balay } 423db4deed7SKarl Rupp } else { 424db4deed7SKarl Rupp for (k=0; k<len; k++) { 4252fa5cd67SKarl Rupp if (x_ptr[k]==0.0) xxt_zero_nnz_0++; 426827bd09bSSatish Balay } 427827bd09bSSatish Balay } 428827bd09bSSatish Balay col_indices[2*i] = off; 429827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = len; 430827bd09bSSatish Balay col_vals[i] = x_ptr; 431827bd09bSSatish Balay } 432db4deed7SKarl Rupp else { 433827bd09bSSatish Balay col_indices[2*i] = 0; 434827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = 0; 435827bd09bSSatish Balay col_vals[i] = x_ptr; 436827bd09bSSatish Balay } 437827bd09bSSatish Balay } 438827bd09bSSatish Balay 439827bd09bSSatish Balay /* close off stages for execution phase */ 440db4deed7SKarl Rupp while (dim!=level) { 441827bd09bSSatish Balay stages[dim++] = i; 44271a0148aSBarry Smith ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 443827bd09bSSatish Balay } 444827bd09bSSatish Balay stages[dim]=i; 445827bd09bSSatish Balay 446827bd09bSSatish Balay xxt_handle->info->n = xxt_handle->mvi->n; 447827bd09bSSatish Balay xxt_handle->info->m = m; 448827bd09bSSatish Balay xxt_handle->info->nnz = xxt_nnz; 449827bd09bSSatish Balay xxt_handle->info->max_nnz = xxt_max_nnz; 450827bd09bSSatish Balay xxt_handle->info->msg_buf_sz = stages[level]-stages[0]; 451a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 452a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar*) malloc(m*sizeof(PetscScalar)); 453827bd09bSSatish Balay xxt_handle->info->x = x; 454827bd09bSSatish Balay xxt_handle->info->col_vals = col_vals; 455827bd09bSSatish Balay xxt_handle->info->col_sz = col_sz; 456827bd09bSSatish Balay xxt_handle->info->col_indices = col_indices; 457827bd09bSSatish Balay xxt_handle->info->stages = stages; 458827bd09bSSatish Balay xxt_handle->info->nsolves = 0; 459827bd09bSSatish Balay xxt_handle->info->tot_solve_time = 0.0; 460827bd09bSSatish Balay 461a501084fSBarry Smith free(segs); 462a501084fSBarry Smith free(u); 463a501084fSBarry Smith free(v); 464a501084fSBarry Smith free(uu); 465a501084fSBarry Smith free(z); 466a501084fSBarry Smith free(w); 467827bd09bSSatish Balay 468827bd09bSSatish Balay return(0); 469827bd09bSSatish Balay } 470827bd09bSSatish Balay 4717b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 4720924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) 473827bd09bSSatish Balay { 47452f87cdaSBarry Smith PetscInt off, len, *iptr; 47552f87cdaSBarry Smith PetscInt level = xxt_handle->level; 47652f87cdaSBarry Smith PetscInt n = xxt_handle->info->n; 47752f87cdaSBarry Smith PetscInt m = xxt_handle->info->m; 47852f87cdaSBarry Smith PetscInt *stages = xxt_handle->info->stages; 47952f87cdaSBarry Smith PetscInt *col_indices = xxt_handle->info->col_indices; 480a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 481a501084fSBarry Smith PetscScalar *solve_uu = xxt_handle->info->solve_uu; 482a501084fSBarry Smith PetscScalar *solve_w = xxt_handle->info->solve_w; 483a501084fSBarry Smith PetscScalar *x = xxt_handle->info->x; 48471a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 485c5df96a5SBarry Smith PetscErrorCode ierr; 486827bd09bSSatish Balay 4870924e98cSBarry Smith PetscFunctionBegin; 488827bd09bSSatish Balay uu_ptr=solve_uu; 489ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 490827bd09bSSatish Balay 491827bd09bSSatish Balay /* x = X.Y^T.b */ 492827bd09bSSatish Balay /* uu = Y^T.b */ 493db4deed7SKarl Rupp for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) { 494c5df96a5SBarry Smith off =*iptr++; 495c5df96a5SBarry Smith len =*iptr++; 496c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 4978b83055fSJed Brown PetscStackCallBLAS("BLASdot",*uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1)); 498827bd09bSSatish Balay } 499827bd09bSSatish Balay 500827bd09bSSatish Balay /* comunication of beta */ 501827bd09bSSatish Balay uu_ptr=solve_uu; 502a00740a5SBarry Smith if (level) {ierr = PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages);CHKERRQ(ierr);} 503827bd09bSSatish Balay 504ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 505827bd09bSSatish Balay 506827bd09bSSatish Balay /* x = X.uu */ 507db4deed7SKarl Rupp for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) { 508c5df96a5SBarry Smith off =*iptr++; 509c5df96a5SBarry Smith len =*iptr++; 510c5df96a5SBarry Smith ierr = PetscBLASIntCast(len,&dlen);CHKERRQ(ierr); 5118b83055fSJed Brown PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1)); 512827bd09bSSatish Balay } 5130924e98cSBarry Smith PetscFunctionReturn(0); 514827bd09bSSatish Balay } 515827bd09bSSatish Balay 5167b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5170924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle) 518827bd09bSSatish Balay { 51952f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 520827bd09bSSatish Balay 5210924e98cSBarry Smith PetscFunctionBegin; 522e32f2f54SBarry Smith if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle); 523827bd09bSSatish Balay 524827bd09bSSatish Balay vals[0]=vals[1]=xxt_handle->id; 525b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 526e32f2f54SBarry 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); 5270924e98cSBarry Smith PetscFunctionReturn(0); 528827bd09bSSatish Balay } 529827bd09bSSatish Balay 5307b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5310924e98cSBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle) 532827bd09bSSatish Balay { 53352f87cdaSBarry Smith PetscInt i, ct, id; 53452f87cdaSBarry Smith PetscInt mask, edge, *iptr; 53552f87cdaSBarry Smith PetscInt *dir, *used; 53652f87cdaSBarry Smith PetscInt sum[4], w[4]; 537a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 53852f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 539a501084fSBarry Smith PetscScalar *lhs, *rhs; 54052f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 541ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle = xxt_handle->mvi->PCTFS_gs_handle; 54252f87cdaSBarry Smith PetscInt *local2global = xxt_handle->mvi->local2global; 54352f87cdaSBarry Smith PetscInt n = xxt_handle->mvi->n; 54452f87cdaSBarry Smith PetscInt m = xxt_handle->mvi->m; 54552f87cdaSBarry Smith PetscInt level = xxt_handle->level; 546ab824b78SBarry Smith PetscInt shared = 0; 547827bd09bSSatish Balay 5480924e98cSBarry Smith PetscFunctionBegin; 54952f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55052f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55152f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55252f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 55352f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 554827bd09bSSatish Balay 555ca8e9878SJed Brown PCTFS_ivec_zero(dir,level+1); 556ca8e9878SJed Brown PCTFS_ivec_zero(nsep,level+1); 557ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 558ca8e9878SJed Brown PCTFS_ivec_set (fo,-1,n+1); 559ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 560827bd09bSSatish Balay 5618cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 5628cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 563827bd09bSSatish Balay 564827bd09bSSatish Balay /* determine the # of unique dof */ 565ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 566ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 567ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 568ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 5693d3eaba7SBarry Smith for (i=0; i<n; i++) { 5702fa5cd67SKarl Rupp if (lhs[i]!=0.0) { 5712fa5cd67SKarl Rupp rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i]; 5722fa5cd67SKarl Rupp } 573827bd09bSSatish Balay } 574b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 575827bd09bSSatish Balay rsum[0]+=0.1; 576827bd09bSSatish Balay rsum[1]+=0.1; 577827bd09bSSatish Balay 57877b4d14cSPeter Brune if (PetscAbsScalar(rsum[0]-rsum[1])>EPS) shared=1; 579827bd09bSSatish Balay 58052f87cdaSBarry Smith xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0]; 58152f87cdaSBarry Smith xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0]; 582827bd09bSSatish Balay 583827bd09bSSatish Balay /* determine separator sets top down */ 5842fa5cd67SKarl Rupp if (shared) { 585db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 586db4deed7SKarl Rupp 587827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 588ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 589ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 590827bd09bSSatish Balay 591827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 592ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 593ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 594827bd09bSSatish Balay 595db4deed7SKarl Rupp for (i=0;i<n;i++) { 596db4deed7SKarl Rupp if (id< mask) { 5972fa5cd67SKarl Rupp if (lhs[i]!=0.0) lhs[i]=1.0; 598827bd09bSSatish Balay } 599db4deed7SKarl Rupp if (id>=mask) { 6002fa5cd67SKarl Rupp if (rhs[i]!=0.0) rhs[i]=1.0; 601827bd09bSSatish Balay } 602827bd09bSSatish Balay } 603827bd09bSSatish Balay 6042fa5cd67SKarl Rupp if (id< mask) PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1); 6052fa5cd67SKarl Rupp else PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1); 606827bd09bSSatish Balay 607827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 608ca8e9878SJed Brown PCTFS_rvec_zero(rsum,4); 609db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 610db4deed7SKarl Rupp if (!used[i]) { 611827bd09bSSatish Balay /* number of unmarked dofs on node */ 612827bd09bSSatish Balay ct++; 613827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 614db4deed7SKarl Rupp if (id< mask) { 615db4deed7SKarl Rupp if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; } 616827bd09bSSatish Balay } 617827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 618db4deed7SKarl Rupp if (id>=mask) { 619db4deed7SKarl Rupp if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; } 620827bd09bSSatish Balay } 621827bd09bSSatish Balay } 622827bd09bSSatish Balay } 623827bd09bSSatish Balay 624827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 625827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 626827bd09bSSatish Balay (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct); 627b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 628b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,4,op,edge); 629827bd09bSSatish Balay rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1; 630827bd09bSSatish Balay 631db4deed7SKarl Rupp if (id<mask) { 632827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 633db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 634db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 635827bd09bSSatish Balay ct++; nfo++; 636827bd09bSSatish Balay 637c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 638827bd09bSSatish Balay 639827bd09bSSatish Balay *--iptr = local2global[i]; 640827bd09bSSatish Balay used[i] = edge; 641827bd09bSSatish Balay } 642827bd09bSSatish Balay } 6432fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 644827bd09bSSatish Balay 645827bd09bSSatish Balay lnsep[edge]=ct; 64652f87cdaSBarry Smith nsep[edge]=(PetscInt) rsum[0]; 647827bd09bSSatish Balay dir [edge]=LEFT; 648827bd09bSSatish Balay } 649827bd09bSSatish Balay 650db4deed7SKarl Rupp if (id>=mask) { 651827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 652db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 653db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 654827bd09bSSatish Balay ct++; nfo++; 655827bd09bSSatish Balay 656c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 657827bd09bSSatish Balay 658827bd09bSSatish Balay *--iptr = local2global[i]; 659827bd09bSSatish Balay used[i] = edge; 660827bd09bSSatish Balay } 661827bd09bSSatish Balay } 6622fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 663827bd09bSSatish Balay 664827bd09bSSatish Balay lnsep[edge] = ct; 66552f87cdaSBarry Smith nsep[edge] = (PetscInt) rsum[1]; 666827bd09bSSatish Balay dir [edge] = RIGHT; 667827bd09bSSatish Balay } 668827bd09bSSatish Balay 669827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 670827bd09bSSatish Balay /* do we need full set of separators for this? */ 671827bd09bSSatish Balay 672827bd09bSSatish Balay /* fold rhs hc into lower */ 6732fa5cd67SKarl Rupp if (id>=mask) id-=mask; 674827bd09bSSatish Balay } 675db4deed7SKarl Rupp } else { 676db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 677827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 678ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 679ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 680827bd09bSSatish Balay 681827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 682ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 683ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 684827bd09bSSatish Balay 685827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 686db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 687db4deed7SKarl Rupp if (!used[i]) { 688827bd09bSSatish Balay /* number of unmarked dofs on node */ 689827bd09bSSatish Balay ct++; 690827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 6912fa5cd67SKarl Rupp if ((id< mask)&&(lhs[i]!=0.0)) sum[0]++; 692827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 6932fa5cd67SKarl Rupp if ((id>=mask)&&(rhs[i]!=0.0)) sum[1]++; 694827bd09bSSatish Balay } 695827bd09bSSatish Balay } 696827bd09bSSatish Balay 697827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 698827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 699b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 700827bd09bSSatish Balay 701827bd09bSSatish Balay /* lhs hc wins */ 702db4deed7SKarl Rupp if (sum[2]>=sum[3]) { 703db4deed7SKarl Rupp if (id<mask) { 704827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 705db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 706db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 707827bd09bSSatish Balay ct++; nfo++; 708827bd09bSSatish Balay *--iptr = local2global[i]; 709827bd09bSSatish Balay used[i]=edge; 710827bd09bSSatish Balay } 711827bd09bSSatish Balay } 7122fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 713827bd09bSSatish Balay lnsep[edge]=ct; 714827bd09bSSatish Balay } 715827bd09bSSatish Balay nsep[edge]=sum[0]; 716827bd09bSSatish Balay dir [edge]=LEFT; 717db4deed7SKarl Rupp } else { /* rhs hc wins */ 718db4deed7SKarl Rupp if (id>=mask) { 719827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 720db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 721db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 722827bd09bSSatish Balay ct++; nfo++; 723827bd09bSSatish Balay *--iptr = local2global[i]; 724827bd09bSSatish Balay used[i]=edge; 725827bd09bSSatish Balay } 726827bd09bSSatish Balay } 7272fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 728827bd09bSSatish Balay lnsep[edge]=ct; 729827bd09bSSatish Balay } 730827bd09bSSatish Balay nsep[edge]=sum[1]; 731827bd09bSSatish Balay dir [edge]=RIGHT; 732827bd09bSSatish Balay } 733827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 734827bd09bSSatish Balay /* do we need full set of separators for this? */ 735827bd09bSSatish Balay 736827bd09bSSatish Balay /* fold rhs hc into lower */ 7372fa5cd67SKarl Rupp if (id>=mask) id-=mask; 738827bd09bSSatish Balay } 739827bd09bSSatish Balay } 740827bd09bSSatish Balay 741827bd09bSSatish Balay 742827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 743db4deed7SKarl Rupp for (ct=i=0; i<n; i++) { 744db4deed7SKarl Rupp if (!used[i]) { 745827bd09bSSatish Balay ct++; nfo++; 746827bd09bSSatish Balay *--iptr = local2global[i]; 747827bd09bSSatish Balay used[i] = edge; 748827bd09bSSatish Balay } 749827bd09bSSatish Balay } 7502fa5cd67SKarl Rupp if (ct>1) PCTFS_ivec_sort(iptr,ct); 751827bd09bSSatish Balay lnsep[edge]=ct; 752827bd09bSSatish Balay nsep [edge]=ct; 753827bd09bSSatish Balay dir [edge]=LEFT; 754827bd09bSSatish Balay 755827bd09bSSatish Balay xxt_handle->info->nsep = nsep; 756827bd09bSSatish Balay xxt_handle->info->lnsep = lnsep; 757827bd09bSSatish Balay xxt_handle->info->fo = fo; 758827bd09bSSatish Balay xxt_handle->info->nfo = nfo; 759827bd09bSSatish Balay 760a501084fSBarry Smith free(dir); 761a501084fSBarry Smith free(lhs); 762a501084fSBarry Smith free(rhs); 763a501084fSBarry Smith free(used); 7640924e98cSBarry Smith PetscFunctionReturn(0); 765827bd09bSSatish Balay } 766827bd09bSSatish Balay 7677b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 7685c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global,PetscInt n,PetscInt m,PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*),void *grid_data) 769827bd09bSSatish Balay { 770827bd09bSSatish Balay mv_info *mvi; 771827bd09bSSatish Balay 772827bd09bSSatish Balay 773a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 774827bd09bSSatish Balay mvi->n = n; 775827bd09bSSatish Balay mvi->m = m; 776827bd09bSSatish Balay mvi->n_global = -1; 777827bd09bSSatish Balay mvi->m_global = -1; 77852f87cdaSBarry Smith mvi->local2global = (PetscInt*)malloc((m+1)*sizeof(PetscInt)); 779ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 780827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 7815c8f6a95SKarl Rupp mvi->matvec = matvec; 782827bd09bSSatish Balay mvi->grid_data = grid_data; 783827bd09bSSatish Balay 784827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 785ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 786827bd09bSSatish Balay 787827bd09bSSatish Balay return(mvi); 788827bd09bSSatish Balay } 789827bd09bSSatish Balay 7907b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 7910924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 792827bd09bSSatish Balay { 7930924e98cSBarry Smith PetscFunctionBegin; 794827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 7950924e98cSBarry Smith PetscFunctionReturn(0); 796827bd09bSSatish Balay } 797827bd09bSSatish Balay 798827bd09bSSatish Balay 799827bd09bSSatish Balay 800