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; 1270e7395fSHong Zhang PetscInt k,start,end,l,row,col,colb,**columnsforrow=coloring->columnsforrow,nz; 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; 1974d3cef9SHong Zhang PetscInt ctype = coloring->ctype,nxloc; 20f5aae955SHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ*)J->data; 21f5aae955SHong Zhang Mat A=aij->A,B=aij->B; 22f5aae955SHong Zhang Mat_SeqAIJ *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data; 23f6af9589SHong Zhang PetscMPIInt rank; 24*573f477fSHong Zhang MatEntry *Jentry=coloring->matentry; 25fcd7ac73SHong Zhang 26fcd7ac73SHong Zhang PetscFunctionBegin; 27f6af9589SHong Zhang ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)J), &rank);CHKERRQ(ierr); 28fcd7ac73SHong Zhang ierr = MatSetUnfactored(J);CHKERRQ(ierr); 29fcd7ac73SHong Zhang ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); 30fcd7ac73SHong Zhang if (flg) { 31fcd7ac73SHong Zhang ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); 32fcd7ac73SHong Zhang } else { 33fcd7ac73SHong Zhang PetscBool assembled; 34fcd7ac73SHong Zhang ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); 35fcd7ac73SHong Zhang if (assembled) { 36fcd7ac73SHong Zhang ierr = MatZeroEntries(J);CHKERRQ(ierr); 37fcd7ac73SHong Zhang } 38fcd7ac73SHong Zhang } 39fcd7ac73SHong Zhang 40f6af9589SHong Zhang if (!coloring->vscale) { //ctype == IS_COLORING_GHOSTED ? 41f6af9589SHong Zhang printf("[%d] create a non-ghostedcoloring->vscale\n",rank); 42f6af9589SHong Zhang ierr = VecDuplicate(x1,&coloring->vscale);CHKERRQ(ierr); 43f6af9589SHong Zhang //SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->scale"); 44fcd7ac73SHong Zhang } 45fcd7ac73SHong Zhang 46fcd7ac73SHong Zhang /* Set w1 = F(x1) */ 47fcd7ac73SHong Zhang if (!coloring->fset) { 48fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 49f6af9589SHong Zhang ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); 50fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 51fcd7ac73SHong Zhang } else { 52fcd7ac73SHong Zhang coloring->fset = PETSC_FALSE; 53fcd7ac73SHong Zhang } 54fcd7ac73SHong Zhang 55fcd7ac73SHong Zhang if (!coloring->w3) { 56f6af9589SHong Zhang ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); 57fcd7ac73SHong Zhang ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); 58fcd7ac73SHong Zhang } 59fcd7ac73SHong Zhang w3 = coloring->w3; 60fcd7ac73SHong Zhang 6174d3cef9SHong Zhang /* Compute vscale = 1./dx - the local scale factors, including ghost points */ 62f6af9589SHong Zhang ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* =J's row ownershipRange, used by ghosted x! */ 63f6af9589SHong Zhang //printf("[%d] start end %d %d\n",rank,start,end); 64f6af9589SHong Zhang 65f6af9589SHong Zhang if (coloring->htype[0] == 'w') { 66f6af9589SHong Zhang /* tacky test; need to make systematic if we add other approaches to computing h*/ 67f6af9589SHong Zhang ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 6874d3cef9SHong Zhang dx = 1.0/((1.0 + unorm)*epsilon); 69f6af9589SHong Zhang ierr = VecSet(coloring->vscale,dx);CHKERRQ(ierr); 7070e7395fSHong Zhang } else { 7174d3cef9SHong Zhang ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); 72f6af9589SHong Zhang ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 73fcd7ac73SHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 7474d3cef9SHong Zhang for (col=0; col<nxloc; col++) { 75fcd7ac73SHong Zhang dx = xx[col]; 7674d3cef9SHong Zhang if (PetscAbsScalar(dx) < umin) { 7774d3cef9SHong Zhang if (PetscRealPart(dx) >= 0.0) dx = umin; 7874d3cef9SHong Zhang else if (PetscRealPart(dx) < 0.0 ) dx = -umin; 7974d3cef9SHong Zhang } 80fcd7ac73SHong Zhang dx *= epsilon; 8174d3cef9SHong Zhang vscale_array[col] = 1.0/dx; 82f6af9589SHong Zhang } 8374d3cef9SHong Zhang ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 84fcd7ac73SHong Zhang ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 8570e7395fSHong Zhang } 86fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 87fcd7ac73SHong Zhang ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88fcd7ac73SHong Zhang ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 89fcd7ac73SHong Zhang } 90fcd7ac73SHong Zhang 911b97d346SHong Zhang /* Loop over each color */ 92f5aae955SHong Zhang nz = 0; 93fcd7ac73SHong Zhang ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 94fcd7ac73SHong Zhang for (k=0; k<coloring->ncolors; k++) { 95fcd7ac73SHong Zhang coloring->currentcolor = k; 96fcd7ac73SHong Zhang 97fcd7ac73SHong Zhang /* 98fcd7ac73SHong Zhang Loop over each column associated with color 99f6af9589SHong Zhang adding the perturbation to the vector w3 = x1 + dx. 100fcd7ac73SHong Zhang */ 101f6af9589SHong Zhang ierr = VecCopy(x1,w3);CHKERRQ(ierr); 102f6af9589SHong Zhang ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); 103f6af9589SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start; 10470e7395fSHong Zhang 105fcd7ac73SHong Zhang for (l=0; l<coloring->ncolumns[k]; l++) { 106f6af9589SHong Zhang col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */ 10770e7395fSHong Zhang /* convert global col to local colb */ 10870e7395fSHong Zhang if (col >= start && col < end) { 10970e7395fSHong Zhang colb = col-start; 11070e7395fSHong Zhang } else { 11170e7395fSHong Zhang #if defined(PETSC_USE_CTABLE) 11270e7395fSHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 11370e7395fSHong Zhang colb--; 11470e7395fSHong Zhang #else 11570e7395fSHong Zhang colb = aij->colmap[col] - 1; /* local column index */ 11670e7395fSHong Zhang #endif 11770e7395fSHong Zhang colb += end - start; 11870e7395fSHong Zhang } 11974d3cef9SHong Zhang w3_array[col] += 1.0/vscale_array[colb]; 120fcd7ac73SHong Zhang } 121fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start; 122fcd7ac73SHong Zhang ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); 123fcd7ac73SHong Zhang 124fcd7ac73SHong Zhang /* 125fcd7ac73SHong Zhang Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) 126fcd7ac73SHong Zhang w2 = F(x1 + dx) - F(x1) 127fcd7ac73SHong Zhang */ 128fcd7ac73SHong Zhang ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 129fcd7ac73SHong Zhang ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); 130fcd7ac73SHong Zhang ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); 131fcd7ac73SHong Zhang ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); 132fcd7ac73SHong Zhang 1331b97d346SHong Zhang /* Loop over rows of vector, putting results into Jacobian matrix */ 134fcd7ac73SHong Zhang ierr = VecGetArray(w2,&y);CHKERRQ(ierr); 135fcd7ac73SHong Zhang for (l=0; l<coloring->nrows[k]; l++) { 136*573f477fSHong Zhang row = Jentry[nz].row; /* local row index */ 137*573f477fSHong Zhang *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col]; 138*573f477fSHong Zhang nz++; 139fcd7ac73SHong Zhang } 140fcd7ac73SHong Zhang ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); 141fcd7ac73SHong Zhang } /* endof for each color */ 142f5aae955SHong 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); 143f5aae955SHong Zhang ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144f5aae955SHong Zhang ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 145fcd7ac73SHong Zhang ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); 146fcd7ac73SHong Zhang 147fcd7ac73SHong Zhang coloring->currentcolor = -1; 148fcd7ac73SHong Zhang PetscFunctionReturn(0); 149fcd7ac73SHong Zhang } 150fcd7ac73SHong Zhang 151fcd7ac73SHong Zhang #undef __FUNCT__ 152fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new" 153fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c) 154fcd7ac73SHong Zhang { 155fcd7ac73SHong Zhang Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 156fcd7ac73SHong Zhang PetscErrorCode ierr; 157fcd7ac73SHong Zhang PetscMPIInt size,*ncolsonproc,*disp,nn; 158fcd7ac73SHong Zhang PetscInt i,n,nrows,j,k,m,ncols,col; 15972c15787SHong Zhang const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL; 160fcd7ac73SHong Zhang PetscInt nis = iscoloring->n,nctot,*cols; 161d3825b63SHong Zhang PetscInt *rowhit,cstart,cend,colb; 162fcd7ac73SHong Zhang IS *isa; 163fcd7ac73SHong Zhang ISLocalToGlobalMapping map = mat->cmap->mapping; 164fcd7ac73SHong Zhang PetscInt ctype=c->ctype; 165fcd7ac73SHong Zhang Mat A=aij->A,B=aij->B; 166fcd7ac73SHong Zhang Mat_SeqAIJ *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data; 167fcd7ac73SHong Zhang PetscScalar *A_val=cspA->a,*B_val=cspB->a; 168fcd7ac73SHong Zhang PetscInt spidx; 1691b97d346SHong Zhang PetscInt *spidxA,*spidxB,nz; 170d3825b63SHong Zhang PetscScalar **valaddrhit; 171*573f477fSHong Zhang MatEntry *Jentry; 172fcd7ac73SHong Zhang 173fcd7ac73SHong Zhang PetscFunctionBegin; 17472c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 17572c15787SHong Zhang if(!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); 17672c15787SHong Zhang ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); 17772c15787SHong Zhang } 178fcd7ac73SHong Zhang 179d3825b63SHong Zhang m = mat->rmap->n; 180fcd7ac73SHong Zhang cstart = mat->cmap->rstart; 181fcd7ac73SHong Zhang cend = mat->cmap->rend; 182fcd7ac73SHong Zhang c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 183fcd7ac73SHong Zhang c->N = mat->cmap->N; 184d3825b63SHong Zhang c->m = m; 185fcd7ac73SHong Zhang c->rstart = mat->rmap->rstart; 186fcd7ac73SHong Zhang 187fcd7ac73SHong Zhang c->ncolors = nis; 188fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 189fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 190fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 191fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 192fcd7ac73SHong Zhang ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 193fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 194fcd7ac73SHong Zhang 195a774a6f1SHong Zhang nz = cspA->nz + cspB->nz; /* total nonzero entries of mat */ 196a774a6f1SHong Zhang ierr = PetscMalloc(nz*sizeof(PetscScalar*),&c->valaddr);CHKERRQ(ierr); 197*573f477fSHong Zhang ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); 198a774a6f1SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(PetscScalar*));CHKERRQ(ierr); 199a774a6f1SHong Zhang 200d3825b63SHong Zhang /* Allow access to data structures of local part of matrix 201d3825b63SHong Zhang - creates aij->colmap which maps global column number to local number in part B */ 202fcd7ac73SHong Zhang if (!aij->colmap) { 203fcd7ac73SHong Zhang ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 204fcd7ac73SHong Zhang } 205fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 206fcd7ac73SHong Zhang ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 207fcd7ac73SHong Zhang 208d3825b63SHong Zhang ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); 209fcd7ac73SHong Zhang nz = 0; 21072c15787SHong Zhang ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 211d3825b63SHong Zhang for (i=0; i<nis; i++) { /* for each local color */ 212fcd7ac73SHong Zhang ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 213fcd7ac73SHong Zhang ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 214fcd7ac73SHong Zhang 215fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 216fcd7ac73SHong Zhang if (n) { 217fcd7ac73SHong Zhang ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 218fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 219fcd7ac73SHong Zhang ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 220fcd7ac73SHong Zhang } else { 221fcd7ac73SHong Zhang c->columns[i] = 0; 222fcd7ac73SHong Zhang } 223fcd7ac73SHong Zhang 224fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 225d3825b63SHong Zhang /* Determine nctot, the total (parallel) number of columns of this color */ 226fcd7ac73SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 227fcd7ac73SHong Zhang ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 228fcd7ac73SHong Zhang 229d3825b63SHong Zhang /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ 230fcd7ac73SHong Zhang ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 231fcd7ac73SHong Zhang ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 232fcd7ac73SHong Zhang nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 233fcd7ac73SHong Zhang if (!nctot) { 234fcd7ac73SHong Zhang ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 235fcd7ac73SHong Zhang } 236fcd7ac73SHong Zhang 237fcd7ac73SHong Zhang disp[0] = 0; 238fcd7ac73SHong Zhang for (j=1; j<size; j++) { 239fcd7ac73SHong Zhang disp[j] = disp[j-1] + ncolsonproc[j-1]; 240fcd7ac73SHong Zhang } 241fcd7ac73SHong Zhang 242d3825b63SHong Zhang /* Get cols, the complete list of columns for this color on each process */ 243fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 244fcd7ac73SHong Zhang ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 245fcd7ac73SHong Zhang ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 246fcd7ac73SHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 247fcd7ac73SHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 248fcd7ac73SHong Zhang nctot = n; 249fcd7ac73SHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 250fcd7ac73SHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 251fcd7ac73SHong Zhang } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 252fcd7ac73SHong Zhang 2531b97d346SHong Zhang /* Mark all rows affect by these columns */ 254d3825b63SHong Zhang ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr); 2551b97d346SHong Zhang for (j=0; j<nctot; j++) { /* loop over columns*/ 256fcd7ac73SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 257fcd7ac73SHong Zhang col = ltog[cols[j]]; 258fcd7ac73SHong Zhang } else { 259fcd7ac73SHong Zhang col = cols[j]; 260fcd7ac73SHong Zhang } 261fcd7ac73SHong Zhang if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */ 262d3825b63SHong Zhang row = A_cj + A_ci[col-cstart]; 263d3825b63SHong Zhang nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; 264fcd7ac73SHong Zhang /* loop over columns of A marking them in rowhit */ 265d3825b63SHong Zhang for (k=0; k<nrows; k++) { 266d3825b63SHong Zhang /* set valaddrhit for part A */ 267fcd7ac73SHong Zhang spidx = spidxA[A_ci[col-cstart] + k]; 268d3825b63SHong Zhang valaddrhit[*row] = &A_val[spidx]; 269a774a6f1SHong Zhang rowhit[*row++] = col - cstart + 1; /* local column index */ 270fcd7ac73SHong Zhang } 271fcd7ac73SHong Zhang } else { /* column is in off-diagonal block of matrix B */ 272fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE) 273fcd7ac73SHong Zhang ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 274fcd7ac73SHong Zhang colb--; 275fcd7ac73SHong Zhang #else 276fcd7ac73SHong Zhang colb = aij->colmap[col] - 1; /* local column index */ 277fcd7ac73SHong Zhang #endif 278fcd7ac73SHong Zhang if (colb == -1) { 279d3825b63SHong Zhang nrows = 0; 280fcd7ac73SHong Zhang } else { 281d3825b63SHong Zhang row = B_cj + B_ci[colb]; 282d3825b63SHong Zhang nrows = B_ci[colb+1] - B_ci[colb]; 283fcd7ac73SHong Zhang } 284fcd7ac73SHong Zhang /* loop over columns of B marking them in rowhit */ 285d3825b63SHong Zhang for (k=0; k<nrows; k++) { 286d3825b63SHong Zhang /* set valaddrhit for part B */ 287fcd7ac73SHong Zhang spidx = spidxB[B_ci[colb] + k]; 288d3825b63SHong Zhang valaddrhit[*row] = &B_val[spidx]; 28970e7395fSHong Zhang rowhit[*row++] = colb + 1 + cend - cstart; /* local column index */ 290fcd7ac73SHong Zhang } 291fcd7ac73SHong Zhang } 292fcd7ac73SHong Zhang } 293fcd7ac73SHong Zhang 294fcd7ac73SHong Zhang /* count the number of hits */ 295fcd7ac73SHong Zhang nrows = 0; 296d3825b63SHong Zhang for (j=0; j<m; j++) { 297fcd7ac73SHong Zhang if (rowhit[j]) nrows++; 298fcd7ac73SHong Zhang } 299fcd7ac73SHong Zhang c->nrows[i] = nrows; 300fcd7ac73SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 301fcd7ac73SHong Zhang ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 302fcd7ac73SHong Zhang ierr = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 303fcd7ac73SHong Zhang nrows = 0; 304d3825b63SHong Zhang for (j=0; j<m; j++) { 305fcd7ac73SHong Zhang if (rowhit[j]) { 306*573f477fSHong Zhang Jentry[nz].row = j; /* local row index */ 307*573f477fSHong Zhang Jentry[nz].col = rowhit[j] - 1; /* local column index */ 308*573f477fSHong Zhang Jentry[nz].valaddr = valaddrhit[j]; /* address of mat value for this entry */ 309*573f477fSHong Zhang 310fcd7ac73SHong Zhang c->rows[i][nrows] = j; /* local row index */ 311a774a6f1SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; /* local column index */ 312*573f477fSHong Zhang c->valaddr[nz] = valaddrhit[j]; /* address of mat value for this entry */ 313*573f477fSHong Zhang nz++; 314fcd7ac73SHong Zhang nrows++; 315fcd7ac73SHong Zhang } 316fcd7ac73SHong Zhang } 317fcd7ac73SHong Zhang ierr = PetscFree(cols);CHKERRQ(ierr); 318fcd7ac73SHong Zhang } 319fcd7ac73SHong 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); 320fcd7ac73SHong Zhang 321a774a6f1SHong Zhang /* create vscale for storing dx */ 322fcd7ac73SHong Zhang if (ctype == IS_COLORING_GLOBAL) { 323d3825b63SHong Zhang ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),m,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 324fcd7ac73SHong Zhang } 325fcd7ac73SHong Zhang 326fcd7ac73SHong Zhang ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 327fcd7ac73SHong Zhang 328*573f477fSHong Zhang c->matentry = Jentry; 329d3825b63SHong Zhang ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr); 330fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); 331fcd7ac73SHong Zhang ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); 33272c15787SHong Zhang if (ctype == IS_COLORING_GHOSTED) { 33372c15787SHong Zhang ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); 33472c15787SHong Zhang } 335fcd7ac73SHong Zhang 336fcd7ac73SHong Zhang mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ; 337fcd7ac73SHong Zhang PetscFunctionReturn(0); 338fcd7ac73SHong Zhang } 339fcd7ac73SHong Zhang 340fcd7ac73SHong Zhang /*------------------------------------------------------*/ 341fcd7ac73SHong Zhang #undef __FUNCT__ 3424a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ" 343dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) 344a64fbb32SBarry Smith { 3456eaac0f3SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 3466849ba73SBarry Smith PetscErrorCode ierr; 347b1d57f15SBarry Smith PetscMPIInt size,*ncolsonproc,*disp,nn; 3481a83f524SJed Brown PetscInt i,n,nrows,j,k,m,ncols,col; 349afcb2eb5SJed Brown const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog; 3501a83f524SJed Brown PetscInt nis = iscoloring->n,nctot,*cols; 351f6d58c54SBarry Smith PetscInt *rowhit,M,cstart,cend,colb; 352b1d57f15SBarry Smith PetscInt *columnsforrow,l; 353b9617806SBarry Smith IS *isa; 354ace3abfcSBarry Smith PetscBool done,flg; 355992144d0SBarry Smith ISLocalToGlobalMapping map = mat->cmap->mapping; 356afcb2eb5SJed Brown PetscInt ctype=c->ctype; 357fcd7ac73SHong Zhang PetscBool new_impl=PETSC_FALSE; 358a64fbb32SBarry Smith 3593a40ed3dSBarry Smith PetscFunctionBegin; 360fcd7ac73SHong Zhang ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr); 361fcd7ac73SHong Zhang if (new_impl){ 362fcd7ac73SHong Zhang ierr = MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr); 363fcd7ac73SHong Zhang PetscFunctionReturn(0); 364fcd7ac73SHong Zhang } 365ce94432eSBarry 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"); 366522c5e43SBarry Smith 367afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr);} 368afcb2eb5SJed Brown else ltog = NULL; 369b9617806SBarry Smith ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 3703acb8795SBarry Smith 371f6d58c54SBarry Smith M = mat->rmap->n; 372f6d58c54SBarry Smith cstart = mat->cmap->rstart; 373f6d58c54SBarry Smith cend = mat->cmap->rend; 374f6d58c54SBarry Smith c->M = mat->rmap->N; /* set the global rows and columns and local rows */ 375f6d58c54SBarry Smith c->N = mat->cmap->N; 376f6d58c54SBarry Smith c->m = mat->rmap->n; 377f6d58c54SBarry Smith c->rstart = mat->rmap->rstart; 378005c665bSBarry Smith 379a64fbb32SBarry Smith c->ncolors = nis; 380b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 381b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); 382b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 383b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr); 384b1d57f15SBarry Smith ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr); 3853bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr); 3866eaac0f3SBarry Smith 3876eaac0f3SBarry Smith /* Allow access to data structures of local part of matrix */ 3886eaac0f3SBarry Smith if (!aij->colmap) { 389ab9863d7SBarry Smith ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); 3906eaac0f3SBarry Smith } 3913acb8795SBarry Smith ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 3923acb8795SBarry Smith ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 3936eaac0f3SBarry Smith 394b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 395b1d57f15SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr); 3966eaac0f3SBarry Smith 397a64fbb32SBarry Smith for (i=0; i<nis; i++) { 398b9b97703SBarry Smith ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 399a64fbb32SBarry Smith ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 4002205254eSKarl Rupp 401fcd7ac73SHong Zhang c->ncolumns[i] = n; /* local number of columns of this color on this process */ 402a64fbb32SBarry Smith if (n) { 403b1d57f15SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); 4043bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); 405b1d57f15SBarry Smith ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); 406a64fbb32SBarry Smith } else { 407a64fbb32SBarry Smith c->columns[i] = 0; 408a64fbb32SBarry Smith } 409a64fbb32SBarry Smith 4108ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 4116eaac0f3SBarry Smith /* Determine the total (parallel) number of columns of this color */ 412ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); 413687f1162SBarry Smith ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); 414b8f8c88eSHong Zhang 4154dc2109aSBarry Smith ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); 416ce94432eSBarry Smith ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4172205254eSKarl Rupp nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j]; 4183a7fca6bSBarry Smith if (!nctot) { 419ae15b995SBarry Smith ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr); 4203a7fca6bSBarry Smith } 4216eaac0f3SBarry Smith 4226eaac0f3SBarry Smith disp[0] = 0; 4236eaac0f3SBarry Smith for (j=1; j<size; j++) { 4246eaac0f3SBarry Smith disp[j] = disp[j-1] + ncolsonproc[j-1]; 4256eaac0f3SBarry Smith } 4266eaac0f3SBarry Smith 4276eaac0f3SBarry Smith /* Get complete list of columns for color on each processor */ 428b1d57f15SBarry Smith ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 429ce94432eSBarry Smith ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 4301d79065fSBarry Smith ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr); 431b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 432b8f8c88eSHong Zhang /* Determine local number of columns of this color on this process, including ghost points */ 433b8f8c88eSHong Zhang nctot = n; 434b8f8c88eSHong Zhang ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr); 435b8f8c88eSHong Zhang ierr = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr); 436f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type"); 4376eaac0f3SBarry Smith 4386eaac0f3SBarry Smith /* 4396eaac0f3SBarry Smith Mark all rows affect by these columns 4406eaac0f3SBarry Smith */ 441b8f8c88eSHong Zhang /* Temporary option to allow for debugging/testing */ 44290d69ab7SBarry Smith flg = PETSC_FALSE; 4430298fd71SBarry Smith ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr); 444f158e583SBarry Smith if (!flg) { /*-----------------------------------------------------------------------------*/ 445f158e583SBarry Smith /* crude, fast version */ 446b1d57f15SBarry Smith ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr); 447a64fbb32SBarry Smith /* loop over columns*/ 4486eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 449b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 450b8f8c88eSHong Zhang col = ltog[cols[j]]; 451b8f8c88eSHong Zhang } else { 4526eaac0f3SBarry Smith col = cols[j]; 453b8f8c88eSHong Zhang } 4546eaac0f3SBarry Smith if (col >= cstart && col < cend) { 4556eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 4566eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 4576eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 4586eaac0f3SBarry Smith } else { 459aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 460cb9801acSJed Brown ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 461fa46199cSSatish Balay colb--; 462b3d2dc96SSatish Balay #else 4636eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 464b3d2dc96SSatish Balay #endif 4656eaac0f3SBarry Smith if (colb == -1) { 4666eaac0f3SBarry Smith m = 0; 4676eaac0f3SBarry Smith } else { 4686eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 4696eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 4706eaac0f3SBarry Smith } 4716eaac0f3SBarry Smith } 472a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 473a64fbb32SBarry Smith for (k=0; k<m; k++) { 474a64fbb32SBarry Smith rowhit[*rows++] = col + 1; 475a64fbb32SBarry Smith } 476a64fbb32SBarry Smith } 4776eaac0f3SBarry Smith 478a64fbb32SBarry Smith /* count the number of hits */ 479a64fbb32SBarry Smith nrows = 0; 4806eaac0f3SBarry Smith for (j=0; j<M; j++) { 481a64fbb32SBarry Smith if (rowhit[j]) nrows++; 482a64fbb32SBarry Smith } 483a64fbb32SBarry Smith c->nrows[i] = nrows; 484b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 485b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 4863bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 487a64fbb32SBarry Smith nrows = 0; 4886eaac0f3SBarry Smith for (j=0; j<M; j++) { 489a64fbb32SBarry Smith if (rowhit[j]) { 490fcd7ac73SHong Zhang c->rows[i][nrows] = j; /* local row index */ 491fcd7ac73SHong Zhang c->columnsforrow[i][nrows] = rowhit[j] - 1; /* global column index */ 492a64fbb32SBarry Smith nrows++; 493a64fbb32SBarry Smith } 494a64fbb32SBarry Smith } 495a64fbb32SBarry Smith } else { /*-------------------------------------------------------------------------------*/ 496f158e583SBarry Smith /* slow version, using rowhit as a linked list */ 497b1d57f15SBarry Smith PetscInt currentcol,fm,mfm; 4986eaac0f3SBarry Smith rowhit[M] = M; 499a64fbb32SBarry Smith nrows = 0; 500a64fbb32SBarry Smith /* loop over columns*/ 5016eaac0f3SBarry Smith for (j=0; j<nctot; j++) { 502b8f8c88eSHong Zhang if (ctype == IS_COLORING_GHOSTED) { 503b8f8c88eSHong Zhang col = ltog[cols[j]]; 504b8f8c88eSHong Zhang } else { 5056eaac0f3SBarry Smith col = cols[j]; 506b8f8c88eSHong Zhang } 5076eaac0f3SBarry Smith if (col >= cstart && col < cend) { 5086eaac0f3SBarry Smith /* column is in diagonal block of matrix */ 5096eaac0f3SBarry Smith rows = A_cj + A_ci[col-cstart]; 5106eaac0f3SBarry Smith m = A_ci[col-cstart+1] - A_ci[col-cstart]; 5116eaac0f3SBarry Smith } else { 512aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 5130f5bd95cSBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 514fa46199cSSatish Balay colb--; 515b3d2dc96SSatish Balay #else 5166eaac0f3SBarry Smith colb = aij->colmap[col] - 1; 517b3d2dc96SSatish Balay #endif 5186eaac0f3SBarry Smith if (colb == -1) { 5196eaac0f3SBarry Smith m = 0; 5206eaac0f3SBarry Smith } else { 5216eaac0f3SBarry Smith rows = B_cj + B_ci[colb]; 5226eaac0f3SBarry Smith m = B_ci[colb+1] - B_ci[colb]; 5236eaac0f3SBarry Smith } 5246eaac0f3SBarry Smith } 525b8f8c88eSHong Zhang 526a64fbb32SBarry Smith /* loop over columns marking them in rowhit */ 5276eaac0f3SBarry Smith fm = M; /* fm points to first entry in linked list */ 528a64fbb32SBarry Smith for (k=0; k<m; k++) { 529a64fbb32SBarry Smith currentcol = *rows++; 530a64fbb32SBarry Smith /* is it already in the list? */ 531a64fbb32SBarry Smith do { 532a64fbb32SBarry Smith mfm = fm; 533a64fbb32SBarry Smith fm = rowhit[fm]; 534a64fbb32SBarry Smith } while (fm < currentcol); 535a64fbb32SBarry Smith /* not in list so add it */ 536a64fbb32SBarry Smith if (fm != currentcol) { 537a64fbb32SBarry Smith nrows++; 538a64fbb32SBarry Smith columnsforrow[currentcol] = col; 539a64fbb32SBarry Smith /* next three lines insert new entry into linked list */ 540a64fbb32SBarry Smith rowhit[mfm] = currentcol; 541a64fbb32SBarry Smith rowhit[currentcol] = fm; 542a64fbb32SBarry Smith fm = currentcol; 543a64fbb32SBarry Smith /* fm points to present position in list since we know the columns are sorted */ 544f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected"); 545a64fbb32SBarry Smith } 546a64fbb32SBarry Smith } 547a64fbb32SBarry Smith c->nrows[i] = nrows; 5482205254eSKarl Rupp 549b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr); 550b1d57f15SBarry Smith ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr); 5513bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr); 552a64fbb32SBarry Smith /* now store the linked list of rows into c->rows[i] */ 553a64fbb32SBarry Smith nrows = 0; 5546eaac0f3SBarry Smith fm = rowhit[M]; 555a64fbb32SBarry Smith do { 556a64fbb32SBarry Smith c->rows[i][nrows] = fm; 557a64fbb32SBarry Smith c->columnsforrow[i][nrows++] = columnsforrow[fm]; 558a64fbb32SBarry Smith fm = rowhit[fm]; 5596eaac0f3SBarry Smith } while (fm < M); 5606eaac0f3SBarry Smith } /* ---------------------------------------------------------------------------------------*/ 561606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 5626eaac0f3SBarry Smith } 56330b34957SBarry Smith 56430b34957SBarry Smith /* Optimize by adding the vscale, and scaleforrow[][] fields */ 56530b34957SBarry Smith /* 56630b34957SBarry Smith vscale will contain the "diagonal" on processor scalings followed by the off processor 56730b34957SBarry Smith */ 5688ee2e534SBarry Smith if (ctype == IS_COLORING_GLOBAL) { 569ce94432eSBarry Smith ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr); 570b1d57f15SBarry Smith ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 57130b34957SBarry Smith for (k=0; k<c->ncolors; k++) { 572b1d57f15SBarry Smith ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 57330b34957SBarry Smith for (l=0; l<c->nrows[k]; l++) { 57430b34957SBarry Smith col = c->columnsforrow[k][l]; 57530b34957SBarry Smith if (col >= cstart && col < cend) { 57630b34957SBarry Smith /* column is in diagonal block of matrix */ 57730b34957SBarry Smith colb = col - cstart; 57830b34957SBarry Smith } else { 57930b34957SBarry Smith /* column is in "off-processor" part */ 58030b34957SBarry Smith #if defined(PETSC_USE_CTABLE) 58130b34957SBarry Smith ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr); 58230b34957SBarry Smith colb--; 58330b34957SBarry Smith #else 58430b34957SBarry Smith colb = aij->colmap[col] - 1; 58530b34957SBarry Smith #endif 58630b34957SBarry Smith colb += cend - cstart; 58730b34957SBarry Smith } 58830b34957SBarry Smith c->vscaleforrow[k][l] = colb; 58930b34957SBarry Smith } 59030b34957SBarry Smith } 591b8f8c88eSHong Zhang } else if (ctype == IS_COLORING_GHOSTED) { 592b8f8c88eSHong Zhang /* Get gtol mapping */ 593afcb2eb5SJed Brown PetscInt N = mat->cmap->N,nlocal,*gtol; 594b8f8c88eSHong Zhang ierr = PetscMalloc((N+1)*sizeof(PetscInt),>ol);CHKERRQ(ierr); 595b8f8c88eSHong Zhang for (i=0; i<N; i++) gtol[i] = -1; 596afcb2eb5SJed Brown ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr); 597afcb2eb5SJed Brown for (i=0; i<nlocal; i++) gtol[ltog[i]] = i; 598b8f8c88eSHong Zhang 599b8f8c88eSHong Zhang c->vscale = 0; /* will be created in MatFDColoringApply() */ 600b8f8c88eSHong Zhang ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr); 601b8f8c88eSHong Zhang for (k=0; k<c->ncolors; k++) { 602b8f8c88eSHong Zhang ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr); 603b8f8c88eSHong Zhang for (l=0; l<c->nrows[k]; l++) { 604b8f8c88eSHong Zhang col = c->columnsforrow[k][l]; /* global column index */ 605b8f8c88eSHong Zhang c->vscaleforrow[k][l] = gtol[col]; /* local column index */ 606b8f8c88eSHong Zhang } 607b8f8c88eSHong Zhang } 608b8f8c88eSHong Zhang ierr = PetscFree(gtol);CHKERRQ(ierr); 609b8f8c88eSHong Zhang } 610b9617806SBarry Smith ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 61130b34957SBarry Smith 612606d414cSSatish Balay ierr = PetscFree(rowhit);CHKERRQ(ierr); 613606d414cSSatish Balay ierr = PetscFree(columnsforrow);CHKERRQ(ierr); 6143acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr); 6153acb8795SBarry Smith ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr); 616afcb2eb5SJed Brown if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr);} 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 618a64fbb32SBarry Smith } 619a64fbb32SBarry Smith 620b9617806SBarry Smith 621b9617806SBarry Smith 622b9617806SBarry Smith 623b9617806SBarry Smith 624b9617806SBarry Smith 625