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); 5952f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle); 6052f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle); 61*5c8f6a95SKarl 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***********************************/ 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 */ 82*5c8f6a95SKarl 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 */ 95*5c8f6a95SKarl 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***********************************/ 1088cda6cd7SBarry Smith PetscInt 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? */ 115db4deed7SKarl Rupp if (b) {PCTFS_rvec_copy(x,b,xxt_handle->mvi->n);} 116827bd09bSSatish Balay do_xxt_solve(xxt_handle,x); 117827bd09bSSatish Balay 118827bd09bSSatish Balay return(0); 119827bd09bSSatish Balay } 120827bd09bSSatish Balay 1217b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 12252f87cdaSBarry Smith PetscInt XXT_free(xxt_ADT xxt_handle) 123827bd09bSSatish Balay { 124827bd09bSSatish Balay 125b1c944f5SJed Brown PCTFS_comm_init(); 126827bd09bSSatish Balay check_handle(xxt_handle); 127827bd09bSSatish Balay n_xxt_handles--; 128827bd09bSSatish Balay 129a501084fSBarry Smith free(xxt_handle->info->nsep); 130a501084fSBarry Smith free(xxt_handle->info->lnsep); 131a501084fSBarry Smith free(xxt_handle->info->fo); 132a501084fSBarry Smith free(xxt_handle->info->stages); 133a501084fSBarry Smith free(xxt_handle->info->solve_uu); 134a501084fSBarry Smith free(xxt_handle->info->solve_w); 135a501084fSBarry Smith free(xxt_handle->info->x); 136a501084fSBarry Smith free(xxt_handle->info->col_vals); 137a501084fSBarry Smith free(xxt_handle->info->col_sz); 138a501084fSBarry Smith free(xxt_handle->info->col_indices); 139a501084fSBarry Smith free(xxt_handle->info); 140a501084fSBarry Smith free(xxt_handle->mvi->local2global); 141ca8e9878SJed Brown PCTFS_gs_free(xxt_handle->mvi->PCTFS_gs_handle); 142a501084fSBarry Smith free(xxt_handle->mvi); 143a501084fSBarry Smith free(xxt_handle); 144827bd09bSSatish Balay 145827bd09bSSatish Balay /* if the check fails we nuke */ 146a501084fSBarry Smith /* if NULL pointer passed to free we nuke */ 147827bd09bSSatish Balay /* if the calls to free fail that's not my problem */ 148827bd09bSSatish Balay return(0); 149827bd09bSSatish Balay } 150827bd09bSSatish Balay 1517b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 15252f87cdaSBarry Smith PetscInt XXT_stats(xxt_ADT xxt_handle) 153827bd09bSSatish Balay { 15452f87cdaSBarry Smith PetscInt op[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD,GL_MIN,GL_MAX,GL_ADD}; 15552f87cdaSBarry Smith PetscInt fop[] = {NON_UNIFORM,GL_MIN,GL_MAX,GL_ADD}; 15652f87cdaSBarry Smith PetscInt vals[9], work[9]; 157a501084fSBarry Smith PetscScalar fvals[3], fwork[3]; 15871a0148aSBarry Smith PetscErrorCode ierr; 159827bd09bSSatish Balay 160b1c944f5SJed Brown PCTFS_comm_init(); 161827bd09bSSatish Balay check_handle(xxt_handle); 162827bd09bSSatish Balay 163827bd09bSSatish Balay /* if factorization not done there are no stats */ 164db4deed7SKarl Rupp if (!xxt_handle->info||!xxt_handle->mvi) { 165b1c944f5SJed Brown if (!PCTFS_my_id) { ierr = PetscPrintf(PETSC_COMM_WORLD,"XXT_stats() :: no stats available!\n");CHKERRQ(ierr); } 166827bd09bSSatish Balay return 1; 167827bd09bSSatish Balay } 168827bd09bSSatish Balay 169827bd09bSSatish Balay vals[0]=vals[1]=vals[2]=xxt_handle->info->nnz; 170827bd09bSSatish Balay vals[3]=vals[4]=vals[5]=xxt_handle->mvi->n; 171827bd09bSSatish Balay vals[6]=vals[7]=vals[8]=xxt_handle->info->msg_buf_sz; 172b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 173827bd09bSSatish Balay 174827bd09bSSatish Balay fvals[0]=fvals[1]=fvals[2] 175827bd09bSSatish Balay =xxt_handle->info->tot_solve_time/xxt_handle->info->nsolves++; 176b1c944f5SJed Brown PCTFS_grop(fvals,fwork,sizeof(fop)/sizeof(fop[0])-1,fop); 177827bd09bSSatish Balay 178db4deed7SKarl Rupp if (!PCTFS_my_id) { 179b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_nnz=%D\n",PCTFS_my_id,vals[0]);CHKERRQ(ierr); 180b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_nnz=%D\n",PCTFS_my_id,vals[1]);CHKERRQ(ierr); 181b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_nnz=%g\n",PCTFS_my_id,1.0*vals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 182b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_nnz=%D\n",PCTFS_my_id,vals[2]);CHKERRQ(ierr); 183b1c944f5SJed 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); 184b1c944f5SJed 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); 185b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_n =%D\n",PCTFS_my_id,vals[3]);CHKERRQ(ierr); 186b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_n =%D\n",PCTFS_my_id,vals[4]);CHKERRQ(ierr); 187b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_n =%g\n",PCTFS_my_id,1.0*vals[5]/PCTFS_num_nodes);CHKERRQ(ierr); 188b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: tot xxt_n =%D\n",PCTFS_my_id,vals[5]);CHKERRQ(ierr); 189b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_buf=%D\n",PCTFS_my_id,vals[6]);CHKERRQ(ierr); 190b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_buf=%D\n",PCTFS_my_id,vals[7]);CHKERRQ(ierr); 191b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_buf=%g\n",PCTFS_my_id,1.0*vals[8]/PCTFS_num_nodes);CHKERRQ(ierr); 192b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: min xxt_slv=%g\n",PCTFS_my_id,fvals[0]);CHKERRQ(ierr); 193b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: max xxt_slv=%g\n",PCTFS_my_id,fvals[1]);CHKERRQ(ierr); 194b1c944f5SJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%D :: avg xxt_slv=%g\n",PCTFS_my_id,fvals[2]/PCTFS_num_nodes);CHKERRQ(ierr); 195827bd09bSSatish Balay } 196827bd09bSSatish Balay 197827bd09bSSatish Balay return(0); 198827bd09bSSatish Balay } 199827bd09bSSatish Balay 200827bd09bSSatish Balay /*************************************xxt.c************************************ 201827bd09bSSatish Balay 202827bd09bSSatish Balay Description: get A_local, local portion of global coarse matrix which 203827bd09bSSatish Balay is a row dist. nxm matrix w/ n<m. 204827bd09bSSatish Balay o my_ml holds address of ML struct associated w/A_local and coarse grid 205827bd09bSSatish Balay o local2global holds global number of column i (i=0,...,m-1) 206827bd09bSSatish Balay o local2global holds global number of row i (i=0,...,n-1) 207827bd09bSSatish Balay o mylocmatvec performs A_local . vec_local (note that gs is performed using 208ca8e9878SJed Brown PCTFS_gs_init/gop). 209827bd09bSSatish Balay 210827bd09bSSatish Balay mylocmatvec = my_ml->Amat[grid_tag].matvec->external; 211827bd09bSSatish Balay mylocmatvec (void :: void *data, double *in, double *out) 212827bd09bSSatish Balay **************************************xxt.c***********************************/ 21352f87cdaSBarry Smith static PetscInt do_xxt_factor(xxt_ADT xxt_handle) 214827bd09bSSatish Balay { 2157b1ae94cSBarry Smith return xxt_generate(xxt_handle); 216827bd09bSSatish Balay } 217827bd09bSSatish Balay 2187b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 21952f87cdaSBarry Smith static PetscInt xxt_generate(xxt_ADT xxt_handle) 220827bd09bSSatish Balay { 22152f87cdaSBarry Smith PetscInt i,j,k,idex; 22252f87cdaSBarry Smith PetscInt dim, col; 223a501084fSBarry Smith PetscScalar *u, *uu, *v, *z, *w, alpha, alpha_w; 22452f87cdaSBarry Smith PetscInt *segs; 22552f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 22652f87cdaSBarry Smith PetscInt off, len; 227a501084fSBarry Smith PetscScalar *x_ptr; 22852f87cdaSBarry Smith PetscInt *iptr, flag; 22952f87cdaSBarry Smith PetscInt start=0, end, work; 23052f87cdaSBarry Smith PetscInt op2[] = {GL_MIN,0}; 231ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle; 23252f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo; 23352f87cdaSBarry Smith PetscInt a_n=xxt_handle->mvi->n; 23452f87cdaSBarry Smith PetscInt a_m=xxt_handle->mvi->m; 23552f87cdaSBarry Smith PetscInt *a_local2global=xxt_handle->mvi->local2global; 23652f87cdaSBarry Smith PetscInt level; 23752f87cdaSBarry Smith PetscInt xxt_nnz=0, xxt_max_nnz=0; 23852f87cdaSBarry Smith PetscInt n, m; 23952f87cdaSBarry Smith PetscInt *col_sz, *col_indices, *stages; 240a501084fSBarry Smith PetscScalar **col_vals, *x; 24152f87cdaSBarry Smith PetscInt n_global; 24252f87cdaSBarry Smith PetscInt xxt_zero_nnz=0; 24352f87cdaSBarry Smith PetscInt xxt_zero_nnz_0=0; 24471a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 245a501084fSBarry Smith PetscScalar dm1 = -1.0; 246d1528f56SBarry Smith PetscErrorCode ierr; 247827bd09bSSatish Balay 248827bd09bSSatish Balay n=xxt_handle->mvi->n; 249827bd09bSSatish Balay nsep=xxt_handle->info->nsep; 250827bd09bSSatish Balay lnsep=xxt_handle->info->lnsep; 251827bd09bSSatish Balay fo=xxt_handle->info->fo; 252827bd09bSSatish Balay end=lnsep[0]; 253827bd09bSSatish Balay level=xxt_handle->level; 254ca8e9878SJed Brown PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle; 255827bd09bSSatish Balay 256827bd09bSSatish Balay /* is there a null space? */ 257827bd09bSSatish Balay /* LATER add in ability to detect null space by checking alpha */ 258827bd09bSSatish Balay for (i=0, j=0; i<=level; i++) 259827bd09bSSatish Balay {j+=nsep[i];} 260827bd09bSSatish Balay 261827bd09bSSatish Balay m = j-xxt_handle->ns; 26222d28d08SBarry Smith if (m!=j) { 26322d28d08SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"xxt_generate() :: null space exists %D %D %D\n",m,j,xxt_handle->ns);CHKERRQ(ierr); 26422d28d08SBarry Smith } 265827bd09bSSatish Balay 266827bd09bSSatish Balay /* get and initialize storage for x local */ 267827bd09bSSatish Balay /* note that x local is nxm and stored by columns */ 26852f87cdaSBarry Smith col_sz = (PetscInt*) malloc(m*sizeof(PetscInt)); 26952f87cdaSBarry Smith col_indices = (PetscInt*) malloc((2*m+1)*sizeof(PetscInt)); 270a501084fSBarry Smith col_vals = (PetscScalar **) malloc(m*sizeof(PetscScalar *)); 271db4deed7SKarl Rupp for (i=j=0; i<m; i++, j+=2) { 272827bd09bSSatish Balay col_indices[j]=col_indices[j+1]=col_sz[i]=-1; 273827bd09bSSatish Balay col_vals[i] = NULL; 274827bd09bSSatish Balay } 275827bd09bSSatish Balay col_indices[j]=-1; 276827bd09bSSatish Balay 277827bd09bSSatish Balay /* size of separators for each sub-hc working from bottom of tree to top */ 278827bd09bSSatish Balay /* this looks like nsep[]=segments */ 27952f87cdaSBarry Smith stages = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 28052f87cdaSBarry Smith segs = (PetscInt*) malloc((level+1)*sizeof(PetscInt)); 281ca8e9878SJed Brown PCTFS_ivec_zero(stages,level+1); 282ca8e9878SJed Brown PCTFS_ivec_copy(segs,nsep,level+1); 283db4deed7SKarl Rupp for (i=0; i<level; i++) { segs[i+1] += segs[i]; } 284827bd09bSSatish Balay stages[0] = segs[0]; 285827bd09bSSatish Balay 286827bd09bSSatish Balay /* temporary vectors */ 287a501084fSBarry Smith u = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 288a501084fSBarry Smith z = (PetscScalar *) malloc(n*sizeof(PetscScalar)); 289a501084fSBarry Smith v = (PetscScalar *) malloc(a_m*sizeof(PetscScalar)); 290a501084fSBarry Smith uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 291a501084fSBarry Smith w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 292827bd09bSSatish Balay 293827bd09bSSatish Balay /* extra nnz due to replication of vertices across separators */ 294db4deed7SKarl Rupp for (i=1, j=0; i<=level; i++) {j+=nsep[i];} 295827bd09bSSatish Balay 296827bd09bSSatish Balay /* storage for sparse x values */ 297827bd09bSSatish Balay n_global = xxt_handle->info->n_global; 298b1c944f5SJed Brown xxt_max_nnz = (PetscInt)(2.5*pow(1.0*n_global,1.6667) + j*n/2)/PCTFS_num_nodes; 299a501084fSBarry Smith x = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar)); 300827bd09bSSatish Balay xxt_nnz = 0; 301827bd09bSSatish Balay 302827bd09bSSatish Balay /* LATER - can embed next sep to fire in gs */ 303827bd09bSSatish Balay /* time to make the donuts - generate X factor */ 304827bd09bSSatish Balay for (dim=i=j=0;i<m;i++) 305827bd09bSSatish Balay { 306827bd09bSSatish Balay /* time to move to the next level? */ 307d4af0d30SBarry Smith while (i==segs[dim]) { 308e32f2f54SBarry Smith if (dim==level) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"dim about to exceed level\n"); 309827bd09bSSatish Balay stages[dim++]=i; 310827bd09bSSatish Balay end+=lnsep[dim]; 311827bd09bSSatish Balay } 312827bd09bSSatish Balay stages[dim]=i; 313827bd09bSSatish Balay 314827bd09bSSatish Balay /* which column are we firing? */ 315827bd09bSSatish Balay /* i.e. set v_l */ 316827bd09bSSatish Balay /* use new seps and do global min across hc to determine which one to fire */ 317827bd09bSSatish Balay (start<end) ? (col=fo[start]) : (col=INT_MAX); 318b1c944f5SJed Brown PCTFS_giop_hc(&col,&work,1,op2,dim); 319827bd09bSSatish Balay 320827bd09bSSatish Balay /* shouldn't need this */ 321db4deed7SKarl Rupp if (col==INT_MAX) { 322f1ed62a8SBarry Smith ierr = PetscInfo(0,"hey ... col==INT_MAX??\n");CHKERRQ(ierr); 323827bd09bSSatish Balay continue; 324827bd09bSSatish Balay } 325827bd09bSSatish Balay 326827bd09bSSatish Balay /* do I own it? I should */ 327ca8e9878SJed Brown PCTFS_rvec_zero(v ,a_m); 328db4deed7SKarl Rupp if (col==fo[start]) { 329827bd09bSSatish Balay start++; 330ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_n); 331e7e72b3dSBarry Smith if (idex!=-1) { 332e7e72b3dSBarry Smith v[idex] = 1.0; j++; 333e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"NOT FOUND!\n"); 334db4deed7SKarl Rupp } else { 335ca8e9878SJed Brown idex=PCTFS_ivec_linear_search(col, a_local2global, a_m); 336827bd09bSSatish Balay if (idex!=-1) 337827bd09bSSatish Balay {v[idex] = 1.0;} 338827bd09bSSatish Balay } 339827bd09bSSatish Balay 340827bd09bSSatish Balay /* perform u = A.v_l */ 341ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 342827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 343827bd09bSSatish Balay 344827bd09bSSatish Balay /* uu = X^T.u_l (local portion) */ 345827bd09bSSatish Balay /* technically only need to zero out first i entries */ 346827bd09bSSatish Balay /* later turn this into an XXT_solve call ? */ 347ca8e9878SJed Brown PCTFS_rvec_zero(uu,m); 348827bd09bSSatish Balay x_ptr=x; 349827bd09bSSatish Balay iptr = col_indices; 350db4deed7SKarl Rupp for (k=0; k<i; k++) { 351827bd09bSSatish Balay off = *iptr++; 352827bd09bSSatish Balay len = *iptr++; 3530805154bSBarry Smith dlen = PetscBLASIntCast(len); 35471a0148aSBarry Smith uu[k] = BLASdot_(&dlen,u+off,&i1,x_ptr,&i1); 355827bd09bSSatish Balay x_ptr+=len; 356827bd09bSSatish Balay } 357827bd09bSSatish Balay 358827bd09bSSatish Balay 359827bd09bSSatish Balay /* uu = X^T.u_l (comm portion) */ 360b1c944f5SJed Brown PCTFS_ssgl_radd (uu, w, dim, stages); 361827bd09bSSatish Balay 362827bd09bSSatish Balay /* z = X.uu */ 363ca8e9878SJed Brown PCTFS_rvec_zero(z,n); 364827bd09bSSatish Balay x_ptr=x; 365827bd09bSSatish Balay iptr = col_indices; 366db4deed7SKarl Rupp for (k=0; k<i; k++) { 367827bd09bSSatish Balay off = *iptr++; 368827bd09bSSatish Balay len = *iptr++; 3690805154bSBarry Smith dlen = PetscBLASIntCast(len); 37071a0148aSBarry Smith BLASaxpy_(&dlen,&uu[k],x_ptr,&i1,z+off,&i1); 371827bd09bSSatish Balay x_ptr+=len; 372827bd09bSSatish Balay } 373827bd09bSSatish Balay 374827bd09bSSatish Balay /* compute v_l = v_l - z */ 375ca8e9878SJed Brown PCTFS_rvec_zero(v+a_n,a_m-a_n); 3760805154bSBarry Smith dlen = PetscBLASIntCast(n); 37771a0148aSBarry Smith BLASaxpy_(&dlen,&dm1,z,&i1,v,&i1); 378827bd09bSSatish Balay 379827bd09bSSatish Balay /* compute u_l = A.v_l */ 380db4deed7SKarl Rupp if (a_n!=a_m) { PCTFS_gs_gop_hc(PCTFS_gs_handle,v,"+\0",dim); } 381ca8e9878SJed Brown PCTFS_rvec_zero(u,n); 382827bd09bSSatish Balay do_matvec(xxt_handle->mvi,v,u); 383827bd09bSSatish Balay 384827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - local portion */ 3850805154bSBarry Smith dlen = PetscBLASIntCast(n); 38671a0148aSBarry Smith alpha = BLASdot_(&dlen,u,&i1,v,&i1); 387827bd09bSSatish Balay /* compute sqrt(alpha) = sqrt(v_l^T.u_l) - comm portion */ 388b1c944f5SJed Brown PCTFS_grop_hc(&alpha, &alpha_w, 1, op, dim); 389827bd09bSSatish Balay 3908f1a2a5eSBarry Smith alpha = (PetscScalar) PetscSqrtReal((PetscReal)alpha); 391827bd09bSSatish Balay 392827bd09bSSatish Balay /* check for small alpha */ 393827bd09bSSatish Balay /* LATER use this to detect and determine null space */ 394c1235816SBarry Smith if (fabs(alpha)<1.0e-14) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"bad alpha! %g\n",alpha); 395827bd09bSSatish Balay 396827bd09bSSatish Balay /* compute v_l = v_l/sqrt(alpha) */ 397ca8e9878SJed Brown PCTFS_rvec_scale(v,1.0/alpha,n); 398827bd09bSSatish Balay 399827bd09bSSatish Balay /* add newly generated column, v_l, to X */ 400827bd09bSSatish Balay flag = 1; 401827bd09bSSatish Balay off=len=0; 402db4deed7SKarl Rupp for (k=0; k<n; k++) { 403db4deed7SKarl Rupp if (v[k]!=0.0) { 404827bd09bSSatish Balay len=k; 405db4deed7SKarl Rupp if (flag) { off=k; flag=0; } 406827bd09bSSatish Balay } 407827bd09bSSatish Balay } 408827bd09bSSatish Balay 409827bd09bSSatish Balay len -= (off-1); 410827bd09bSSatish Balay 411db4deed7SKarl Rupp if (len>0) { 412db4deed7SKarl Rupp if ((xxt_nnz+len)>xxt_max_nnz) { 413f1ed62a8SBarry Smith ierr = PetscInfo(0,"increasing space for X by 2x!\n");CHKERRQ(ierr); 414827bd09bSSatish Balay xxt_max_nnz *= 2; 415a501084fSBarry Smith x_ptr = (PetscScalar *) malloc(xxt_max_nnz*sizeof(PetscScalar)); 416ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,x,xxt_nnz); 417a501084fSBarry Smith free(x); 418827bd09bSSatish Balay x = x_ptr; 419827bd09bSSatish Balay x_ptr+=xxt_nnz; 420827bd09bSSatish Balay } 421827bd09bSSatish Balay xxt_nnz += len; 422ca8e9878SJed Brown PCTFS_rvec_copy(x_ptr,v+off,len); 423827bd09bSSatish Balay 424827bd09bSSatish Balay /* keep track of number of zeros */ 425db4deed7SKarl Rupp if (dim) { 426db4deed7SKarl Rupp for (k=0; k<len; k++) { 427db4deed7SKarl Rupp if (x_ptr[k]==0.0) { xxt_zero_nnz++; } 428827bd09bSSatish Balay } 429db4deed7SKarl Rupp } else { 430db4deed7SKarl Rupp for (k=0; k<len; k++) { 431db4deed7SKarl Rupp if (x_ptr[k]==0.0) {xxt_zero_nnz_0++;} 432827bd09bSSatish Balay } 433827bd09bSSatish Balay } 434827bd09bSSatish Balay col_indices[2*i] = off; 435827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = len; 436827bd09bSSatish Balay col_vals[i] = x_ptr; 437827bd09bSSatish Balay } 438db4deed7SKarl Rupp else { 439827bd09bSSatish Balay col_indices[2*i] = 0; 440827bd09bSSatish Balay col_sz[i] = col_indices[2*i+1] = 0; 441827bd09bSSatish Balay col_vals[i] = x_ptr; 442827bd09bSSatish Balay } 443827bd09bSSatish Balay } 444827bd09bSSatish Balay 445827bd09bSSatish Balay /* close off stages for execution phase */ 446db4deed7SKarl Rupp while (dim!=level) { 447827bd09bSSatish Balay stages[dim++]=i; 44871a0148aSBarry Smith ierr = PetscInfo2(0,"disconnected!!! dim(%D)!=level(%D)\n",dim,level);CHKERRQ(ierr); 449827bd09bSSatish Balay } 450827bd09bSSatish Balay stages[dim]=i; 451827bd09bSSatish Balay 452827bd09bSSatish Balay xxt_handle->info->n=xxt_handle->mvi->n; 453827bd09bSSatish Balay xxt_handle->info->m=m; 454827bd09bSSatish Balay xxt_handle->info->nnz=xxt_nnz; 455827bd09bSSatish Balay xxt_handle->info->max_nnz=xxt_max_nnz; 456827bd09bSSatish Balay xxt_handle->info->msg_buf_sz=stages[level]-stages[0]; 457a501084fSBarry Smith xxt_handle->info->solve_uu = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 458a501084fSBarry Smith xxt_handle->info->solve_w = (PetscScalar *) malloc(m*sizeof(PetscScalar)); 459827bd09bSSatish Balay xxt_handle->info->x=x; 460827bd09bSSatish Balay xxt_handle->info->col_vals=col_vals; 461827bd09bSSatish Balay xxt_handle->info->col_sz=col_sz; 462827bd09bSSatish Balay xxt_handle->info->col_indices=col_indices; 463827bd09bSSatish Balay xxt_handle->info->stages=stages; 464827bd09bSSatish Balay xxt_handle->info->nsolves=0; 465827bd09bSSatish Balay xxt_handle->info->tot_solve_time=0.0; 466827bd09bSSatish Balay 467a501084fSBarry Smith free(segs); 468a501084fSBarry Smith free(u); 469a501084fSBarry Smith free(v); 470a501084fSBarry Smith free(uu); 471a501084fSBarry Smith free(z); 472a501084fSBarry Smith free(w); 473827bd09bSSatish Balay 474827bd09bSSatish Balay return(0); 475827bd09bSSatish Balay } 476827bd09bSSatish Balay 4777b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 4780924e98cSBarry Smith static PetscErrorCode do_xxt_solve(xxt_ADT xxt_handle, PetscScalar *uc) 479827bd09bSSatish Balay { 48052f87cdaSBarry Smith PetscInt off, len, *iptr; 48152f87cdaSBarry Smith PetscInt level =xxt_handle->level; 48252f87cdaSBarry Smith PetscInt n =xxt_handle->info->n; 48352f87cdaSBarry Smith PetscInt m =xxt_handle->info->m; 48452f87cdaSBarry Smith PetscInt *stages =xxt_handle->info->stages; 48552f87cdaSBarry Smith PetscInt *col_indices=xxt_handle->info->col_indices; 486a501084fSBarry Smith PetscScalar *x_ptr, *uu_ptr; 487a501084fSBarry Smith PetscScalar *solve_uu=xxt_handle->info->solve_uu; 488a501084fSBarry Smith PetscScalar *solve_w =xxt_handle->info->solve_w; 489a501084fSBarry Smith PetscScalar *x =xxt_handle->info->x; 49071a0148aSBarry Smith PetscBLASInt i1 = 1,dlen; 491827bd09bSSatish Balay 4920924e98cSBarry Smith PetscFunctionBegin; 493827bd09bSSatish Balay uu_ptr=solve_uu; 494ca8e9878SJed Brown PCTFS_rvec_zero(uu_ptr,m); 495827bd09bSSatish Balay 496827bd09bSSatish Balay /* x = X.Y^T.b */ 497827bd09bSSatish Balay /* uu = Y^T.b */ 498db4deed7SKarl Rupp for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) { 499827bd09bSSatish Balay off=*iptr++; len=*iptr++; 5000805154bSBarry Smith dlen = PetscBLASIntCast(len); 50171a0148aSBarry Smith *uu_ptr++ = BLASdot_(&dlen,uc+off,&i1,x_ptr,&i1); 502827bd09bSSatish Balay } 503827bd09bSSatish Balay 504827bd09bSSatish Balay /* comunication of beta */ 505827bd09bSSatish Balay uu_ptr=solve_uu; 506b1c944f5SJed Brown if (level) { PCTFS_ssgl_radd(uu_ptr, solve_w, level, stages); } 507827bd09bSSatish Balay 508ca8e9878SJed Brown PCTFS_rvec_zero(uc,n); 509827bd09bSSatish Balay 510827bd09bSSatish Balay /* x = X.uu */ 511db4deed7SKarl Rupp for (x_ptr=x,iptr=col_indices; *iptr!=-1; x_ptr+=len) { 512827bd09bSSatish Balay off=*iptr++; len=*iptr++; 5130805154bSBarry Smith dlen = PetscBLASIntCast(len); 51471a0148aSBarry Smith BLASaxpy_(&dlen,uu_ptr++,x_ptr,&i1,uc+off,&i1); 515827bd09bSSatish Balay } 5160924e98cSBarry Smith PetscFunctionReturn(0); 517827bd09bSSatish Balay } 518827bd09bSSatish Balay 5197b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5200924e98cSBarry Smith static PetscErrorCode check_handle(xxt_ADT xxt_handle) 521827bd09bSSatish Balay { 52252f87cdaSBarry Smith PetscInt vals[2], work[2], op[] = {NON_UNIFORM,GL_MIN,GL_MAX}; 523827bd09bSSatish Balay 5240924e98cSBarry Smith PetscFunctionBegin; 525e32f2f54SBarry Smith if (!xxt_handle) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"check_handle() :: bad handle :: NULL %D\n",xxt_handle); 526827bd09bSSatish Balay 527827bd09bSSatish Balay vals[0]=vals[1]=xxt_handle->id; 528b1c944f5SJed Brown PCTFS_giop(vals,work,sizeof(op)/sizeof(op[0])-1,op); 529e32f2f54SBarry 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); 5300924e98cSBarry Smith PetscFunctionReturn(0); 531827bd09bSSatish Balay } 532827bd09bSSatish Balay 5337b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 5340924e98cSBarry Smith static PetscErrorCode det_separators(xxt_ADT xxt_handle) 535827bd09bSSatish Balay { 53652f87cdaSBarry Smith PetscInt i, ct, id; 53752f87cdaSBarry Smith PetscInt mask, edge, *iptr; 53852f87cdaSBarry Smith PetscInt *dir, *used; 53952f87cdaSBarry Smith PetscInt sum[4], w[4]; 540a501084fSBarry Smith PetscScalar rsum[4], rw[4]; 54152f87cdaSBarry Smith PetscInt op[] = {GL_ADD,0}; 542a501084fSBarry Smith PetscScalar *lhs, *rhs; 54352f87cdaSBarry Smith PetscInt *nsep, *lnsep, *fo, nfo=0; 544ca8e9878SJed Brown PCTFS_gs_ADT PCTFS_gs_handle=xxt_handle->mvi->PCTFS_gs_handle; 54552f87cdaSBarry Smith PetscInt *local2global=xxt_handle->mvi->local2global; 54652f87cdaSBarry Smith PetscInt n=xxt_handle->mvi->n; 54752f87cdaSBarry Smith PetscInt m=xxt_handle->mvi->m; 54852f87cdaSBarry Smith PetscInt level=xxt_handle->level; 549ab824b78SBarry Smith PetscInt shared=0; 550827bd09bSSatish Balay 5510924e98cSBarry Smith PetscFunctionBegin; 55252f87cdaSBarry Smith dir = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55352f87cdaSBarry Smith nsep = (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55452f87cdaSBarry Smith lnsep= (PetscInt*)malloc(sizeof(PetscInt)*(level+1)); 55552f87cdaSBarry Smith fo = (PetscInt*)malloc(sizeof(PetscInt)*(n+1)); 55652f87cdaSBarry Smith used = (PetscInt*)malloc(sizeof(PetscInt)*n); 557827bd09bSSatish Balay 558ca8e9878SJed Brown PCTFS_ivec_zero(dir ,level+1); 559ca8e9878SJed Brown PCTFS_ivec_zero(nsep ,level+1); 560ca8e9878SJed Brown PCTFS_ivec_zero(lnsep,level+1); 561ca8e9878SJed Brown PCTFS_ivec_set (fo ,-1,n+1); 562ca8e9878SJed Brown PCTFS_ivec_zero(used,n); 563827bd09bSSatish Balay 5648cda6cd7SBarry Smith lhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 5658cda6cd7SBarry Smith rhs = (PetscScalar*)malloc(sizeof(PetscScalar)*m); 566827bd09bSSatish Balay 567827bd09bSSatish Balay /* determine the # of unique dof */ 568ca8e9878SJed Brown PCTFS_rvec_zero(lhs,m); 569ca8e9878SJed Brown PCTFS_rvec_set(lhs,1.0,n); 570ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",level); 571ca8e9878SJed Brown PCTFS_rvec_zero(rsum,2); 572db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 573827bd09bSSatish Balay if (lhs[i]!=0.0) 574827bd09bSSatish Balay {rsum[0]+=1.0/lhs[i]; rsum[1]+=lhs[i];} 575827bd09bSSatish Balay } 576b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,2,op,level); 577827bd09bSSatish Balay rsum[0]+=0.1; 578827bd09bSSatish Balay rsum[1]+=0.1; 579827bd09bSSatish Balay 580db4deed7SKarl Rupp if (fabs(rsum[0]-rsum[1])>EPS) { shared=1; } 581827bd09bSSatish Balay 58252f87cdaSBarry Smith xxt_handle->info->n_global=xxt_handle->info->m_global=(PetscInt) rsum[0]; 58352f87cdaSBarry Smith xxt_handle->mvi->n_global =xxt_handle->mvi->m_global =(PetscInt) rsum[0]; 584827bd09bSSatish Balay 585827bd09bSSatish Balay /* determine separator sets top down */ 586827bd09bSSatish Balay if (shared) 587827bd09bSSatish Balay { 588db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 589db4deed7SKarl Rupp 590827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 591ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 592ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 593827bd09bSSatish Balay 594827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 595ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 596ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 597827bd09bSSatish Balay 598db4deed7SKarl Rupp for (i=0;i<n;i++) { 599db4deed7SKarl Rupp if (id< mask) { 600db4deed7SKarl Rupp if (lhs[i]!=0.0) { lhs[i]=1.0; } 601827bd09bSSatish Balay } 602db4deed7SKarl Rupp 603db4deed7SKarl Rupp if (id>=mask) { 604db4deed7SKarl Rupp if (rhs[i]!=0.0) { rhs[i]=1.0; } 605827bd09bSSatish Balay } 606827bd09bSSatish Balay } 607827bd09bSSatish Balay 608db4deed7SKarl Rupp if (id< mask) { PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge-1); } 609db4deed7SKarl Rupp else { PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge-1); } 610827bd09bSSatish Balay 611827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 612ca8e9878SJed Brown PCTFS_rvec_zero(rsum,4); 613db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 614db4deed7SKarl Rupp 615db4deed7SKarl Rupp if (!used[i]) { 616827bd09bSSatish Balay /* number of unmarked dofs on node */ 617827bd09bSSatish Balay ct++; 618827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 619db4deed7SKarl Rupp if (id< mask) { 620db4deed7SKarl Rupp if (lhs[i]!=0.0) { sum[0]++; rsum[0]+=1.0/lhs[i]; } 621827bd09bSSatish Balay } 622827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 623db4deed7SKarl Rupp if (id>=mask) { 624db4deed7SKarl Rupp if (rhs[i]!=0.0) { sum[1]++; rsum[1]+=1.0/rhs[i]; } 625827bd09bSSatish Balay } 626827bd09bSSatish Balay } 627db4deed7SKarl Rupp 628827bd09bSSatish Balay } 629827bd09bSSatish Balay 630827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 631827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 632827bd09bSSatish Balay (id<mask) ? (rsum[2]=ct) : (rsum[3]=ct); 633b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 634b1c944f5SJed Brown PCTFS_grop_hc(rsum,rw,4,op,edge); 635827bd09bSSatish Balay rsum[0]+=0.1; rsum[1]+=0.1; rsum[2]+=0.1; rsum[3]+=0.1; 636827bd09bSSatish Balay 637db4deed7SKarl Rupp if (id<mask) { 638827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 639db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 640db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 641827bd09bSSatish Balay ct++; nfo++; 642827bd09bSSatish Balay 643c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 644827bd09bSSatish Balay 645827bd09bSSatish Balay *--iptr = local2global[i]; 646827bd09bSSatish Balay used[i]=edge; 647827bd09bSSatish Balay } 648827bd09bSSatish Balay } 649ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 650827bd09bSSatish Balay 651827bd09bSSatish Balay lnsep[edge]=ct; 65252f87cdaSBarry Smith nsep[edge]=(PetscInt) rsum[0]; 653827bd09bSSatish Balay dir [edge]=LEFT; 654827bd09bSSatish Balay } 655827bd09bSSatish Balay 656db4deed7SKarl Rupp if (id>=mask) { 657827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 658db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 659db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 660827bd09bSSatish Balay ct++; nfo++; 661827bd09bSSatish Balay 662c1235816SBarry Smith if (nfo>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"nfo about to exceed n\n"); 663827bd09bSSatish Balay 664827bd09bSSatish Balay *--iptr = local2global[i]; 665827bd09bSSatish Balay used[i]=edge; 666827bd09bSSatish Balay } 667827bd09bSSatish Balay } 668ca8e9878SJed Brown if (ct>1) { PCTFS_ivec_sort(iptr,ct); } 669827bd09bSSatish Balay 670827bd09bSSatish Balay lnsep[edge]=ct; 67152f87cdaSBarry Smith nsep[edge]= (PetscInt) rsum[1]; 672827bd09bSSatish Balay dir [edge]=RIGHT; 673827bd09bSSatish Balay } 674827bd09bSSatish Balay 675827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 676827bd09bSSatish Balay /* do we need full set of separators for this? */ 677827bd09bSSatish Balay 678827bd09bSSatish Balay /* fold rhs hc into lower */ 679db4deed7SKarl Rupp if (id>=mask) { id-=mask; } 680827bd09bSSatish Balay } 681db4deed7SKarl Rupp } else { 682db4deed7SKarl Rupp for (iptr=fo+n,id=PCTFS_my_id,mask=PCTFS_num_nodes>>1,edge=level;edge>0;edge--,mask>>=1) { 683827bd09bSSatish Balay /* set rsh of hc, fire, and collect lhs responses */ 684ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_zero(lhs,m) : PCTFS_rvec_set(lhs,1.0,m); 685ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,lhs,"+\0",edge); 686827bd09bSSatish Balay 687827bd09bSSatish Balay /* set lsh of hc, fire, and collect rhs responses */ 688ca8e9878SJed Brown (id<mask) ? PCTFS_rvec_set(rhs,1.0,m) : PCTFS_rvec_zero(rhs,m); 689ca8e9878SJed Brown PCTFS_gs_gop_hc(PCTFS_gs_handle,rhs,"+\0",edge); 690827bd09bSSatish Balay 691827bd09bSSatish Balay /* count number of dofs I own that have signal and not in sep set */ 692db4deed7SKarl Rupp for (PCTFS_ivec_zero(sum,4),ct=i=0;i<n;i++) { 693db4deed7SKarl Rupp if (!used[i]) { 694827bd09bSSatish Balay /* number of unmarked dofs on node */ 695827bd09bSSatish Balay ct++; 696827bd09bSSatish Balay /* number of dofs to be marked on lhs hc */ 697827bd09bSSatish Balay if ((id< mask)&&(lhs[i]!=0.0)) { sum[0]++; } 698827bd09bSSatish Balay /* number of dofs to be marked on rhs hc */ 699827bd09bSSatish Balay if ((id>=mask)&&(rhs[i]!=0.0)) { sum[1]++; } 700827bd09bSSatish Balay } 701827bd09bSSatish Balay } 702827bd09bSSatish Balay 703827bd09bSSatish Balay /* go for load balance - choose half with most unmarked dofs, bias LHS */ 704827bd09bSSatish Balay (id<mask) ? (sum[2]=ct) : (sum[3]=ct); 705b1c944f5SJed Brown PCTFS_giop_hc(sum,w,4,op,edge); 706827bd09bSSatish Balay 707827bd09bSSatish Balay /* lhs hc wins */ 708db4deed7SKarl Rupp if (sum[2]>=sum[3]) { 709db4deed7SKarl Rupp if (id<mask) { 710827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 711db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 712db4deed7SKarl Rupp if ((!used[i])&&(lhs[i]!=0.0)) { 713827bd09bSSatish Balay ct++; nfo++; 714827bd09bSSatish Balay *--iptr = local2global[i]; 715827bd09bSSatish Balay used[i]=edge; 716827bd09bSSatish Balay } 717827bd09bSSatish Balay } 718ca8e9878SJed Brown if (ct>1) { PCTFS_ivec_sort(iptr,ct); } 719827bd09bSSatish Balay lnsep[edge]=ct; 720827bd09bSSatish Balay } 721827bd09bSSatish Balay nsep[edge]=sum[0]; 722827bd09bSSatish Balay dir [edge]=LEFT; 723db4deed7SKarl Rupp } else { /* rhs hc wins */ 724db4deed7SKarl Rupp if (id>=mask) { 725827bd09bSSatish Balay /* mark dofs I own that have signal and not in sep set */ 726db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 727db4deed7SKarl Rupp if ((!used[i])&&(rhs[i]!=0.0)) { 728827bd09bSSatish Balay ct++; nfo++; 729827bd09bSSatish Balay *--iptr = local2global[i]; 730827bd09bSSatish Balay used[i]=edge; 731827bd09bSSatish Balay } 732827bd09bSSatish Balay } 733ca8e9878SJed Brown if (ct>1) {PCTFS_ivec_sort(iptr,ct);} 734827bd09bSSatish Balay lnsep[edge]=ct; 735827bd09bSSatish Balay } 736827bd09bSSatish Balay nsep[edge]=sum[1]; 737827bd09bSSatish Balay dir [edge]=RIGHT; 738827bd09bSSatish Balay } 739827bd09bSSatish Balay /* LATER or we can recur on these to order seps at this level */ 740827bd09bSSatish Balay /* do we need full set of separators for this? */ 741827bd09bSSatish Balay 742827bd09bSSatish Balay /* fold rhs hc into lower */ 743db4deed7SKarl Rupp if (id>=mask) { id-=mask; } 744827bd09bSSatish Balay } 745827bd09bSSatish Balay } 746827bd09bSSatish Balay 747827bd09bSSatish Balay 748827bd09bSSatish Balay /* level 0 is on processor case - so mark the remainder */ 749db4deed7SKarl Rupp for (ct=i=0;i<n;i++) { 750db4deed7SKarl Rupp if (!used[i]) { 751827bd09bSSatish Balay ct++; nfo++; 752827bd09bSSatish Balay *--iptr = local2global[i]; 753827bd09bSSatish Balay used[i]=edge; 754827bd09bSSatish Balay } 755827bd09bSSatish Balay } 756ca8e9878SJed Brown if (ct>1) { PCTFS_ivec_sort(iptr,ct); } 757827bd09bSSatish Balay lnsep[edge]=ct; 758827bd09bSSatish Balay nsep [edge]=ct; 759827bd09bSSatish Balay dir [edge]=LEFT; 760827bd09bSSatish Balay 761827bd09bSSatish Balay xxt_handle->info->nsep=nsep; 762827bd09bSSatish Balay xxt_handle->info->lnsep=lnsep; 763827bd09bSSatish Balay xxt_handle->info->fo=fo; 764827bd09bSSatish Balay xxt_handle->info->nfo=nfo; 765827bd09bSSatish Balay 766a501084fSBarry Smith free(dir); 767a501084fSBarry Smith free(lhs); 768a501084fSBarry Smith free(rhs); 769a501084fSBarry Smith free(used); 7700924e98cSBarry Smith PetscFunctionReturn(0); 771827bd09bSSatish Balay } 772827bd09bSSatish Balay 7737b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 774*5c8f6a95SKarl Rupp static mv_info *set_mvi(PetscInt *local2global,PetscInt n,PetscInt m,PetscErrorCode (*matvec)(mv_info*,PetscScalar*,PetscScalar*),void *grid_data) 775827bd09bSSatish Balay { 776827bd09bSSatish Balay mv_info *mvi; 777827bd09bSSatish Balay 778827bd09bSSatish Balay 779a501084fSBarry Smith mvi = (mv_info*)malloc(sizeof(mv_info)); 780827bd09bSSatish Balay mvi->n=n; 781827bd09bSSatish Balay mvi->m=m; 782827bd09bSSatish Balay mvi->n_global=-1; 783827bd09bSSatish Balay mvi->m_global=-1; 78452f87cdaSBarry Smith mvi->local2global=(PetscInt*)malloc((m+1)*sizeof(PetscInt)); 785ca8e9878SJed Brown PCTFS_ivec_copy(mvi->local2global,local2global,m); 786827bd09bSSatish Balay mvi->local2global[m] = INT_MAX; 787*5c8f6a95SKarl Rupp mvi->matvec=matvec; 788827bd09bSSatish Balay mvi->grid_data=grid_data; 789827bd09bSSatish Balay 790827bd09bSSatish Balay /* set xxt communication handle to perform restricted matvec */ 791ca8e9878SJed Brown mvi->PCTFS_gs_handle = PCTFS_gs_init(local2global, m, PCTFS_num_nodes); 792827bd09bSSatish Balay 793827bd09bSSatish Balay return(mvi); 794827bd09bSSatish Balay } 795827bd09bSSatish Balay 7967b1ae94cSBarry Smith /**************************************xxt.c***********************************/ 7970924e98cSBarry Smith static PetscErrorCode do_matvec(mv_info *A, PetscScalar *v, PetscScalar *u) 798827bd09bSSatish Balay { 7990924e98cSBarry Smith PetscFunctionBegin; 800827bd09bSSatish Balay A->matvec((mv_info*)A->grid_data,v,u); 8010924e98cSBarry Smith PetscFunctionReturn(0); 802827bd09bSSatish Balay } 803827bd09bSSatish Balay 804827bd09bSSatish Balay 805827bd09bSSatish Balay 806