1a64fbb32SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 3a64fbb32SBarry Smith 4ab9863d7SBarry Smith extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 5a64fbb32SBarry Smith 64a2ae208SSatish Balay #undef __FUNCT__ 7fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringApply_MPIAIJ" 8fcd7ac73SHong Zhang PetscErrorCode MatFDColoringApply_MPIAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) 9fcd7ac73SHong Zhang { 10fcd7ac73SHong Zhang PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; 11fcd7ac73SHong Zhang PetscErrorCode ierr; 12*f5aae955SHong Zhang PetscInt k,start,end,l,row,col,**vscaleforrow; 13fcd7ac73SHong Zhang PetscScalar dx,*y,*xx,*w3_array; 14fcd7ac73SHong Zhang PetscScalar *vscale_array; 15fcd7ac73SHong Zhang PetscReal epsilon = coloring->error_rel,umin = coloring->umin,unorm; 16fcd7ac73SHong Zhang Vec w1 = coloring->w1,w2=coloring->w2,w3; 17fcd7ac73SHong Zhang void *fctx = coloring->fctx; 18fcd7ac73SHong Zhang PetscBool flg = PETSC_FALSE; 19fcd7ac73SHong Zhang PetscInt ctype = coloring->ctype,N,col_start=0,col_end=0; 20fcd7ac73SHong Zhang Vec x1_tmp; 21*f5aae955SHong Zhang PetscInt nz; 22*f5aae955SHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ*)J->data; 23*f5aae955SHong Zhang Mat A=aij->A,B=aij->B; 24*f5aae955SHong Zhang Mat_SeqAIJ *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data; 25fcd7ac73SHong Zhang 26fcd7ac73SHong Zhang PetscFunctionBegin; 27fcd7ac73SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 28fcd7ac73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 29fcd7ac73SHong Zhang if (flg) { 30fcd7ac73SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 31fcd7ac73SHong Zhang } else { 32fcd7ac73SHong Zhang PetscBool assembled; 33fcd7ac73SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 34fcd7ac73SHong Zhang if (assembled) { 35fcd7ac73SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 36fcd7ac73SHong Zhang } 37fcd7ac73SHong Zhang } 38fcd7ac73SHong Zhang 39fcd7ac73SHong Zhang x1_tmp = x1; 40fcd7ac73SHong Zhang if (!coloring->vscale) { 41fcd7ac73SHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr); 42fcd7ac73SHong Zhang } 43fcd7ac73SHong Zhang 44fcd7ac73SHong Zhang if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/ 45fcd7ac73SHong Zhang ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr); 46fcd7ac73SHong Zhang } 47fcd7ac73SHong Zhang ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */ 48fcd7ac73SHong Zhang 49fcd7ac73SHong Zhang /* Set w1 = F(x1) */ 50fcd7ac73SHong Zhang if (!coloring->fset) { 51fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 52fcd7ac73SHong Zhang ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr); 53fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 54fcd7ac73SHong Zhang } else { 55fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 56fcd7ac73SHong Zhang } 57fcd7ac73SHong Zhang 58fcd7ac73SHong Zhang if (!coloring->w3) { 59fcd7ac73SHong Zhang ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr); 60fcd7ac73SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 61fcd7ac73SHong Zhang } 62fcd7ac73SHong Zhang w3 = coloring->w3; 63fcd7ac73SHong Zhang 64fcd7ac73SHong Zhang /* Compute all the local scale factors, including ghost points */ 65fcd7ac73SHong Zhang ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr); 66fcd7ac73SHong Zhang ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr); 67fcd7ac73SHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 68fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 69fcd7ac73SHong Zhang col_start = 0; col_end = N; 70fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GLOBAL) { 71fcd7ac73SHong Zhang xx = xx - start; 72fcd7ac73SHong Zhang vscale_array = vscale_array - start; 73fcd7ac73SHong Zhang col_start = start; col_end = N + start; 74fcd7ac73SHong Zhang } 75fcd7ac73SHong Zhang for (col=col_start; col<col_end; col++) { 76fcd7ac73SHong Zhang /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */ 77fcd7ac73SHong Zhang if (coloring->htype[0] == 'w') { 78fcd7ac73SHong Zhang dx = 1.0 + unorm; 79fcd7ac73SHong Zhang } else { 80fcd7ac73SHong Zhang dx = xx[col]; 81fcd7ac73SHong Zhang } 82fcd7ac73SHong Zhang if (dx == (PetscScalar)0.0) dx = 1.0; 83fcd7ac73SHong Zhang if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 84fcd7ac73SHong Zhang else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 85fcd7ac73SHong Zhang dx *= epsilon; 86fcd7ac73SHong Zhang vscale_array[col] = (PetscScalar)1.0/dx; 87fcd7ac73SHong Zhang } 88fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start; 89fcd7ac73SHong Zhang ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 90fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 91fcd7ac73SHong Zhang ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 92fcd7ac73SHong Zhang ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 93fcd7ac73SHong Zhang } 94fcd7ac73SHong Zhang 95fcd7ac73SHong Zhang if (coloring->vscaleforrow) { 96fcd7ac73SHong Zhang vscaleforrow = coloring->vscaleforrow; 97fcd7ac73SHong Zhang } else SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow"); 98fcd7ac73SHong Zhang 99fcd7ac73SHong Zhang /* 100fcd7ac73SHong Zhang Loop over each color 101fcd7ac73SHong Zhang */ 102*f5aae955SHong Zhang nz = 0; 103fcd7ac73SHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 104fcd7ac73SHong Zhang for (k=0; k<coloring->ncolors; k++) { 105fcd7ac73SHong Zhang coloring->currentcolor = k; 106fcd7ac73SHong Zhang 107fcd7ac73SHong Zhang ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr); 108fcd7ac73SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 109fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 110fcd7ac73SHong Zhang /* 111fcd7ac73SHong Zhang Loop over each column associated with color 112fcd7ac73SHong Zhang adding the perturbation to the vector w3. 113fcd7ac73SHong Zhang */ 114fcd7ac73SHong Zhang for (l=0; l<coloring->ncolumns[k]; l++) { 115fcd7ac73SHong Zhang col = coloring->columns[k][l]; /* local column of the matrix we are probing for */ 116fcd7ac73SHong Zhang if (coloring->htype[0] == 'w') { 117fcd7ac73SHong Zhang dx = 1.0 + unorm; 118fcd7ac73SHong Zhang } else { 119fcd7ac73SHong Zhang dx = xx[col]; 120fcd7ac73SHong Zhang } 121fcd7ac73SHong Zhang if (dx == (PetscScalar)0.0) dx = 1.0; 122fcd7ac73SHong Zhang if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; 123fcd7ac73SHong Zhang else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; 124fcd7ac73SHong Zhang dx *= epsilon; 125fcd7ac73SHong Zhang if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter"); 126fcd7ac73SHong Zhang w3_array[col] += dx; 127fcd7ac73SHong Zhang } 128fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 129fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 130fcd7ac73SHong Zhang 131fcd7ac73SHong Zhang /* 132fcd7ac73SHong Zhang Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 133fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 134fcd7ac73SHong Zhang */ 135fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 136fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 137fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 138fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 139fcd7ac73SHong Zhang 140fcd7ac73SHong Zhang /* 141fcd7ac73SHong Zhang Loop over rows of vector, putting results into Jacobian matrix 142fcd7ac73SHong Zhang */ 143fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 144fcd7ac73SHong Zhang for (l=0; l<coloring->nrows[k]; l++) { 145fcd7ac73SHong Zhang row = coloring->rows[k][l]; /* local row index */ 146*f5aae955SHong Zhang *(coloring->spaddr[nz++]) = vscale_array[vscaleforrow[k][l]]*y[row]; 147fcd7ac73SHong Zhang } 148fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 149fcd7ac73SHong Zhang } /* endof for each color */ 150*f5aae955SHong Zhang if (nz != cspA->nz + cspB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,cspA->nz+cspB->nz); 151*f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152*f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153fcd7ac73SHong Zhang 154fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) xx = xx + start; 155fcd7ac73SHong Zhang ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 156fcd7ac73SHong Zhang ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr); 157fcd7ac73SHong Zhang 158fcd7ac73SHong Zhang coloring->currentcolor = -1; 159fcd7ac73SHong Zhang PetscFunctionReturn(0); 160fcd7ac73SHong Zhang } 161fcd7ac73SHong Zhang 162fcd7ac73SHong Zhang #undef __FUNCT__ 163fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new" 164fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c) 165fcd7ac73SHong Zhang { 166fcd7ac73SHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 167fcd7ac73SHong Zhang PetscErrorCode ierr; 168fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 169fcd7ac73SHong Zhang PetscInt i,n,nrows,j,k,m,ncols,col; 170fcd7ac73SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog; 171fcd7ac73SHong Zhang PetscInt nis = iscoloring->n,nctot,*cols; 172fcd7ac73SHong Zhang PetscInt *rowhit,M,cstart,cend,colb; 173fcd7ac73SHong Zhang PetscInt *columnsforrow,l; 174fcd7ac73SHong Zhang IS *isa; 175fcd7ac73SHong Zhang ISLocalToGlobalMapping map = mat->cmap->mapping; 176fcd7ac73SHong Zhang PetscInt ctype=c->ctype; 177fcd7ac73SHong Zhang Mat A=aij->A,B=aij->B; 178fcd7ac73SHong Zhang Mat_SeqAIJ *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data; 179fcd7ac73SHong Zhang PetscScalar *A_val=cspA->a,*B_val=cspB->a; 180fcd7ac73SHong Zhang PetscInt spidx; 181fcd7ac73SHong Zhang 182fcd7ac73SHong Zhang PetscFunctionBegin; 183fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 184fcd7ac73SHong Zhang 185fcd7ac73SHong Zhang if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr);} 186fcd7ac73SHong Zhang else ltog = NULL; 187fcd7ac73SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 188fcd7ac73SHong Zhang 189fcd7ac73SHong Zhang M = mat->rmap->n; 190fcd7ac73SHong Zhang cstart = mat->cmap->rstart; 191fcd7ac73SHong Zhang cend = mat->cmap->rend; 192fcd7ac73SHong Zhang c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 193fcd7ac73SHong Zhang c->N = mat->cmap->N; 194fcd7ac73SHong Zhang c->m = mat->rmap->n; 195fcd7ac73SHong Zhang c->rstart = mat->rmap->rstart; 196fcd7ac73SHong Zhang 197fcd7ac73SHong Zhang c->ncolors = nis; 198fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 199fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 200fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 201fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 202fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 203fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 204fcd7ac73SHong Zhang 205fcd7ac73SHong Zhang /* Allow access to data structures of local part of matrix */ 206fcd7ac73SHong Zhang if (!aij->colmap) { 207fcd7ac73SHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 208fcd7ac73SHong Zhang } 209fcd7ac73SHong Zhang PetscInt *spidxA,*spidxB,nz; 210fcd7ac73SHong Zhang PetscScalar **spaddrhit; 211fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 212fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 213fcd7ac73SHong Zhang 214fcd7ac73SHong Zhang ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 215fcd7ac73SHong Zhang ierr = PetscMalloc((M+1)*sizeof(PetscScalar*),&spaddrhit);CHKERRQ(ierr); 216fcd7ac73SHong Zhang ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 217fcd7ac73SHong Zhang nz = cspA->nz + cspB->nz; /* total nonzero entries of mat */ 218fcd7ac73SHong Zhang ierr = PetscMalloc(nz*sizeof(PetscScalar*),&c->spaddr);CHKERRQ(ierr); 219fcd7ac73SHong Zhang 220fcd7ac73SHong Zhang nz = 0; 221fcd7ac73SHong Zhang for (i=0; i<nis; i++) { 222fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 223fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 224fcd7ac73SHong Zhang 225fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 226fcd7ac73SHong Zhang if (n) { 227fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 228fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 229fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 230fcd7ac73SHong Zhang } else { 231fcd7ac73SHong Zhang c->columns[i] = 0; 232fcd7ac73SHong Zhang } 233fcd7ac73SHong Zhang 234fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 235fcd7ac73SHong Zhang /* Determine the total (parallel) number of columns of this color */ 236fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 237fcd7ac73SHong Zhang ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 238fcd7ac73SHong Zhang 239fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 240fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 241fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 242fcd7ac73SHong Zhang if (!nctot) { 243fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 244fcd7ac73SHong Zhang } 245fcd7ac73SHong Zhang 246fcd7ac73SHong Zhang disp[0] = 0; 247fcd7ac73SHong Zhang for (j=1; j<size; j++) { 248fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 249fcd7ac73SHong Zhang } 250fcd7ac73SHong Zhang 251fcd7ac73SHong Zhang /* Get complete list of columns for color on each processor */ 252fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 253fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 254fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 255fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 256fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 257fcd7ac73SHong Zhang nctot = n; 258fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 259fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 260fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 261fcd7ac73SHong Zhang 262fcd7ac73SHong Zhang /* 263fcd7ac73SHong Zhang Mark all rows affect by these columns 264fcd7ac73SHong Zhang */ 265fcd7ac73SHong Zhang /*-----------------------------------------------------------------------------*/ 266fcd7ac73SHong Zhang /* crude, fast version */ 267fcd7ac73SHong Zhang ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 268fcd7ac73SHong Zhang /* loop over columns*/ 269fcd7ac73SHong Zhang for (j=0; j<nctot; j++) { 270fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 271fcd7ac73SHong Zhang col = ltog[cols[j]]; 272fcd7ac73SHong Zhang } else { 273fcd7ac73SHong Zhang col = cols[j]; 274fcd7ac73SHong Zhang } 275fcd7ac73SHong Zhang if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */ 276fcd7ac73SHong Zhang rows = A_cj + A_ci[col-cstart]; 277fcd7ac73SHong Zhang m = A_ci[col-cstart+1] - A_ci[col-cstart]; 278fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 279fcd7ac73SHong Zhang for (k=0; k<m; k++) { 280fcd7ac73SHong Zhang /* set spaddrhit for part A */ 281fcd7ac73SHong Zhang spidx = spidxA[A_ci[col-cstart] + k]; 282fcd7ac73SHong Zhang spaddrhit[*rows] = &A_val[spidx]; 283fcd7ac73SHong Zhang rowhit[*rows++] = col + 1; 284fcd7ac73SHong Zhang } 285fcd7ac73SHong Zhang } else { /* column is in off-diagonal block of matrix B */ 286fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 287fcd7ac73SHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 288fcd7ac73SHong Zhang colb--; 289fcd7ac73SHong Zhang #else 290fcd7ac73SHong Zhang colb = aij->colmap[col] - 1; /* local column index */ 291fcd7ac73SHong Zhang #endif 292fcd7ac73SHong Zhang if (colb == -1) { 293fcd7ac73SHong Zhang m = 0; 294fcd7ac73SHong Zhang } else { 295fcd7ac73SHong Zhang rows = B_cj + B_ci[colb]; 296fcd7ac73SHong Zhang m = B_ci[colb+1] - B_ci[colb]; 297fcd7ac73SHong Zhang } 298fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 299fcd7ac73SHong Zhang for (k=0; k<m; k++) { 300fcd7ac73SHong Zhang /* set spaddrhit for part B */ 301fcd7ac73SHong Zhang spidx = spidxB[B_ci[colb] + k]; 302fcd7ac73SHong Zhang spaddrhit[*rows] = &B_val[spidx]; 303fcd7ac73SHong Zhang rowhit[*rows++] = col + 1; 304fcd7ac73SHong Zhang } 305fcd7ac73SHong Zhang } 306fcd7ac73SHong Zhang } 307fcd7ac73SHong Zhang 308fcd7ac73SHong Zhang /* count the number of hits */ 309fcd7ac73SHong Zhang nrows = 0; 310fcd7ac73SHong Zhang for (j=0; j<M; j++) { 311fcd7ac73SHong Zhang if (rowhit[j]) nrows++; 312fcd7ac73SHong Zhang } 313fcd7ac73SHong Zhang c->nrows[i] = nrows; 314fcd7ac73SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 315fcd7ac73SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 316fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 317fcd7ac73SHong Zhang nrows = 0; 318fcd7ac73SHong Zhang for (j=0; j<M; j++) { 319fcd7ac73SHong Zhang if (rowhit[j]) { 320fcd7ac73SHong Zhang c->rows[i][nrows] = j; /* local row index */ 321fcd7ac73SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; /* global column index */ 322fcd7ac73SHong Zhang c->spaddr[nz++] = spaddrhit[j]; /* address of mat value for this entry */ 323fcd7ac73SHong Zhang nrows++; 324fcd7ac73SHong Zhang } 325fcd7ac73SHong Zhang } 326fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 327fcd7ac73SHong Zhang } 328fcd7ac73SHong Zhang if (nz != cspA->nz + cspB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,cspA->nz+cspB->nz); 329fcd7ac73SHong Zhang 330fcd7ac73SHong Zhang /* Optimize by adding the vscale, and scaleforrow[][] fields */ 331fcd7ac73SHong Zhang /* 332fcd7ac73SHong Zhang vscale will contain the "diagonal" on processor scalings followed by the off processor 333fcd7ac73SHong Zhang */ 334fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 335fcd7ac73SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 336fcd7ac73SHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 337fcd7ac73SHong Zhang for (k=0; k<c->ncolors; k++) { 338fcd7ac73SHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 339fcd7ac73SHong Zhang for (l=0; l<c->nrows[k]; l++) { 340fcd7ac73SHong Zhang col = c->columnsforrow[k][l]; 341fcd7ac73SHong Zhang if (col >= cstart && col < cend) { 342fcd7ac73SHong Zhang /* column is in diagonal block of matrix */ 343fcd7ac73SHong Zhang colb = col - cstart; 344fcd7ac73SHong Zhang } else { 345fcd7ac73SHong Zhang /* column is in "off-processor" part */ 346fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 347fcd7ac73SHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 348fcd7ac73SHong Zhang colb--; 349fcd7ac73SHong Zhang #else 350fcd7ac73SHong Zhang colb = aij->colmap[col] - 1; 351fcd7ac73SHong Zhang #endif 352fcd7ac73SHong Zhang colb += cend - cstart; 353fcd7ac73SHong Zhang } 354fcd7ac73SHong Zhang c->vscaleforrow[k][l] = colb; 355fcd7ac73SHong Zhang } 356fcd7ac73SHong Zhang } 357fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 358fcd7ac73SHong Zhang /* Get gtol mapping */ 359fcd7ac73SHong Zhang PetscInt N = mat->cmap->N,nlocal,*gtol; 360fcd7ac73SHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 361fcd7ac73SHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 362fcd7ac73SHong Zhang ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr); 363fcd7ac73SHong Zhang for (i=0; i<nlocal; i++) gtol[ltog[i]] = i; 364fcd7ac73SHong Zhang 365fcd7ac73SHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 366fcd7ac73SHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 367fcd7ac73SHong Zhang for (k=0; k<c->ncolors; k++) { 368fcd7ac73SHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 369fcd7ac73SHong Zhang for (l=0; l<c->nrows[k]; l++) { 370fcd7ac73SHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 371fcd7ac73SHong Zhang 372fcd7ac73SHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 373fcd7ac73SHong Zhang } 374fcd7ac73SHong Zhang } 375fcd7ac73SHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 376fcd7ac73SHong Zhang } 377fcd7ac73SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 378fcd7ac73SHong Zhang 379fcd7ac73SHong Zhang ierr = PetscFree(rowhit);CHKERRQ(ierr); 380fcd7ac73SHong Zhang ierr = PetscFree(spaddrhit);CHKERRQ(ierr); 381fcd7ac73SHong Zhang ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 382fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 383fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 384fcd7ac73SHong Zhang if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr);} 385fcd7ac73SHong Zhang 386fcd7ac73SHong Zhang mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ; 387fcd7ac73SHong Zhang PetscFunctionReturn(0); 388fcd7ac73SHong Zhang } 389fcd7ac73SHong Zhang 390fcd7ac73SHong Zhang /*------------------------------------------------------*/ 391fcd7ac73SHong Zhang #undef __FUNCT__ 3924a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ" 393dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 394a64fbb32SBarry Smith { 3956eaac0f3SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 3966849ba73SBarry Smith PetscErrorCode ierr; 397b1d57f15SBarry Smith PetscMPIInt size,*ncolsonproc,*disp,nn; 3981a83f524SJed Brown PetscInt i,n,nrows,j,k,m,ncols,col; 399afcb2eb5SJed Brown const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog; 4001a83f524SJed Brown PetscInt nis = iscoloring->n,nctot,*cols; 401f6d58c54SBarry Smith PetscInt *rowhit,M,cstart,cend,colb; 402b1d57f15SBarry Smith PetscInt *columnsforrow,l; 403b9617806SBarry Smith IS *isa; 404ace3abfcSBarry Smith PetscBool done,flg; 405992144d0SBarry Smith ISLocalToGlobalMapping map = mat->cmap->mapping; 406afcb2eb5SJed Brown PetscInt ctype=c->ctype; 407fcd7ac73SHong Zhang PetscBool new_impl=PETSC_FALSE; 408a64fbb32SBarry Smith 4093a40ed3dSBarry Smith PetscFunctionBegin; 410fcd7ac73SHong Zhang ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr); 411fcd7ac73SHong Zhang if (new_impl){ 412fcd7ac73SHong Zhang ierr = MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr); 413fcd7ac73SHong Zhang PetscFunctionReturn(0); 414fcd7ac73SHong Zhang } 415ce94432eSBarry Smith if (ctype == IS_COLORING_GHOSTED && !map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 416522c5e43SBarry Smith 417afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr);} 418afcb2eb5SJed Brown else ltog = NULL; 419b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 4203acb8795SBarry Smith 421f6d58c54SBarry Smith M = mat->rmap->n; 422f6d58c54SBarry Smith cstart = mat->cmap->rstart; 423f6d58c54SBarry Smith cend = mat->cmap->rend; 424f6d58c54SBarry Smith c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 425f6d58c54SBarry Smith c->N = mat->cmap->N; 426f6d58c54SBarry Smith c->m = mat->rmap->n; 427f6d58c54SBarry Smith c->rstart = mat->rmap->rstart; 428005c665bSBarry Smith 429a64fbb32SBarry Smith c->ncolors = nis; 430b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 431b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 432b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 433b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 434b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 4353bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 4366eaac0f3SBarry Smith 4376eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 4386eaac0f3SBarry Smith if (!aij->colmap) { 439ab9863d7SBarry Smith ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 4406eaac0f3SBarry Smith } 4413acb8795SBarry Smith ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 4423acb8795SBarry Smith ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 4436eaac0f3SBarry Smith 444b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 445b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 4466eaac0f3SBarry Smith 447a64fbb32SBarry Smith for (i=0; i<nis; i++) { 448b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 449a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 4502205254eSKarl Rupp 451fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 452a64fbb32SBarry Smith if (n) { 453b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 4543bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 455b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 456a64fbb32SBarry Smith } else { 457a64fbb32SBarry Smith c->columns[i] = 0; 458a64fbb32SBarry Smith } 459a64fbb32SBarry Smith 4608ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 4616eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 462ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 463687f1162SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 464b8f8c88eSHong Zhang 4654dc2109aSBarry Smith ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 466ce94432eSBarry Smith ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4672205254eSKarl Rupp nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 4683a7fca6bSBarry Smith if (!nctot) { 469ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 4703a7fca6bSBarry Smith } 4716eaac0f3SBarry Smith 4726eaac0f3SBarry Smith disp[0] = 0; 4736eaac0f3SBarry Smith for (j=1; j<size; j++) { 4746eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 4756eaac0f3SBarry Smith } 4766eaac0f3SBarry Smith 4776eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 478b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 479ce94432eSBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4801d79065fSBarry Smith ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 481b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 482b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 483b8f8c88eSHong Zhang nctot = n; 484b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 485b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 486f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 4876eaac0f3SBarry Smith 4886eaac0f3SBarry Smith /* 4896eaac0f3SBarry Smith Mark all rows affect by these columns 4906eaac0f3SBarry Smith */ 491b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 49290d69ab7SBarry Smith flg = PETSC_FALSE; 4930298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr); 494f158e583SBarry Smith if (!flg) { /*-----------------------------------------------------------------------------*/ 495f158e583SBarry Smith /* crude, fast version */ 496b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 497a64fbb32SBarry Smith /* loop over columns*/ 4986eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 499b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 500b8f8c88eSHong Zhang col = ltog[cols[j]]; 501b8f8c88eSHong Zhang } else { 5026eaac0f3SBarry Smith col = cols[j]; 503b8f8c88eSHong Zhang } 5046eaac0f3SBarry Smith if (col >= cstart && col < cend) { 5056eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 5066eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 5076eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 5086eaac0f3SBarry Smith } else { 509aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 510cb9801acSJed Brown ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 511fa46199cSSatish Balay colb--; 512b3d2dc96SSatish Balay #else 5136eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 514b3d2dc96SSatish Balay #endif 5156eaac0f3SBarry Smith if (colb == -1) { 5166eaac0f3SBarry Smith m = 0; 5176eaac0f3SBarry Smith } else { 5186eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 5196eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 5206eaac0f3SBarry Smith } 5216eaac0f3SBarry Smith } 522a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 523a64fbb32SBarry Smith for (k=0; k<m; k++) { 524a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 525a64fbb32SBarry Smith } 526a64fbb32SBarry Smith } 5276eaac0f3SBarry Smith 528a64fbb32SBarry Smith /* count the number of hits */ 529a64fbb32SBarry Smith nrows = 0; 5306eaac0f3SBarry Smith for (j=0; j<M; j++) { 531a64fbb32SBarry Smith if (rowhit[j]) nrows++; 532a64fbb32SBarry Smith } 533a64fbb32SBarry Smith c->nrows[i] = nrows; 534b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 535b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 5363bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 537a64fbb32SBarry Smith nrows = 0; 5386eaac0f3SBarry Smith for (j=0; j<M; j++) { 539a64fbb32SBarry Smith if (rowhit[j]) { 540fcd7ac73SHong Zhang c->rows[i][nrows] = j; /* local row index */ 541fcd7ac73SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; /* global column index */ 542a64fbb32SBarry Smith nrows++; 543a64fbb32SBarry Smith } 544a64fbb32SBarry Smith } 545a64fbb32SBarry Smith } else { /*-------------------------------------------------------------------------------*/ 546f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 547b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 5486eaac0f3SBarry Smith rowhit[M] = M; 549a64fbb32SBarry Smith nrows = 0; 550a64fbb32SBarry Smith /* loop over columns*/ 5516eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 552b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 553b8f8c88eSHong Zhang col = ltog[cols[j]]; 554b8f8c88eSHong Zhang } else { 5556eaac0f3SBarry Smith col = cols[j]; 556b8f8c88eSHong Zhang } 5576eaac0f3SBarry Smith if (col >= cstart && col < cend) { 5586eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 5596eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 5606eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 5616eaac0f3SBarry Smith } else { 562aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 5630f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 564fa46199cSSatish Balay colb--; 565b3d2dc96SSatish Balay #else 5666eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 567b3d2dc96SSatish Balay #endif 5686eaac0f3SBarry Smith if (colb == -1) { 5696eaac0f3SBarry Smith m = 0; 5706eaac0f3SBarry Smith } else { 5716eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 5726eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 5736eaac0f3SBarry Smith } 5746eaac0f3SBarry Smith } 575b8f8c88eSHong Zhang 576a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 5776eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 578a64fbb32SBarry Smith for (k=0; k<m; k++) { 579a64fbb32SBarry Smith currentcol = *rows++; 580a64fbb32SBarry Smith /* is it already in the list? */ 581a64fbb32SBarry Smith do { 582a64fbb32SBarry Smith mfm = fm; 583a64fbb32SBarry Smith fm = rowhit[fm]; 584a64fbb32SBarry Smith } while (fm < currentcol); 585a64fbb32SBarry Smith /* not in list so add it */ 586a64fbb32SBarry Smith if (fm != currentcol) { 587a64fbb32SBarry Smith nrows++; 588a64fbb32SBarry Smith columnsforrow[currentcol] = col; 589a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 590a64fbb32SBarry Smith rowhit[mfm] = currentcol; 591a64fbb32SBarry Smith rowhit[currentcol] = fm; 592a64fbb32SBarry Smith fm = currentcol; 593a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 594f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 595a64fbb32SBarry Smith } 596a64fbb32SBarry Smith } 597a64fbb32SBarry Smith c->nrows[i] = nrows; 5982205254eSKarl Rupp 599b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 600b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 6013bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 602a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 603a64fbb32SBarry Smith nrows = 0; 6046eaac0f3SBarry Smith fm = rowhit[M]; 605a64fbb32SBarry Smith do { 606a64fbb32SBarry Smith c->rows[i][nrows] = fm; 607a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 608a64fbb32SBarry Smith fm = rowhit[fm]; 6096eaac0f3SBarry Smith } while (fm < M); 6106eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 611606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 6126eaac0f3SBarry Smith } 61330b34957SBarry Smith 61430b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 61530b34957SBarry Smith /* 61630b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 61730b34957SBarry Smith */ 6188ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 619ce94432eSBarry Smith ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 620b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 62130b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 622b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 62330b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 62430b34957SBarry Smith col = c->columnsforrow[k][l]; 62530b34957SBarry Smith if (col >= cstart && col < cend) { 62630b34957SBarry Smith /* column is in diagonal block of matrix */ 62730b34957SBarry Smith colb = col - cstart; 62830b34957SBarry Smith } else { 62930b34957SBarry Smith /* column is in "off-processor" part */ 63030b34957SBarry Smith #if defined(PETSC_USE_CTABLE) 63130b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 63230b34957SBarry Smith colb--; 63330b34957SBarry Smith #else 63430b34957SBarry Smith colb = aij->colmap[col] - 1; 63530b34957SBarry Smith #endif 63630b34957SBarry Smith colb += cend - cstart; 63730b34957SBarry Smith } 63830b34957SBarry Smith c->vscaleforrow[k][l] = colb; 63930b34957SBarry Smith } 64030b34957SBarry Smith } 641b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 642b8f8c88eSHong Zhang /* Get gtol mapping */ 643afcb2eb5SJed Brown PetscInt N = mat->cmap->N,nlocal,*gtol; 644b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 645b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 646afcb2eb5SJed Brown ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr); 647afcb2eb5SJed Brown for (i=0; i<nlocal; i++) gtol[ltog[i]] = i; 648b8f8c88eSHong Zhang 649b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 650b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 651b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 652b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 653b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 654b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 6552205254eSKarl Rupp 656b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 657b8f8c88eSHong Zhang } 658b8f8c88eSHong Zhang } 659b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 660b8f8c88eSHong Zhang } 661b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 66230b34957SBarry Smith 663606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 664606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 6653acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 6663acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 667afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr);} 6683a40ed3dSBarry Smith PetscFunctionReturn(0); 669a64fbb32SBarry Smith } 670a64fbb32SBarry Smith 671b9617806SBarry Smith 672b9617806SBarry Smith 673b9617806SBarry Smith 674b9617806SBarry Smith 675b9617806SBarry Smith 676