xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision f6af95899a7fa5afca16675068e593b249b86d93)
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*f6af9589SHong Zhang   PetscInt       k,start,end,l,row,col,**vscaleforrow,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;
19fcd7ac73SHong Zhang   PetscInt       ctype   = coloring->ctype,N,col_start=0,col_end=0;
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;
23*f6af9589SHong Zhang   PetscMPIInt    rank;
24fcd7ac73SHong Zhang 
25fcd7ac73SHong Zhang   PetscFunctionBegin;
26*f6af9589SHong Zhang   //ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);
27*f6af9589SHong Zhang   //ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD);
28*f6af9589SHong Zhang   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)J), &rank);CHKERRQ(ierr);
29fcd7ac73SHong Zhang   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
30fcd7ac73SHong Zhang   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
31fcd7ac73SHong Zhang   if (flg) {
32fcd7ac73SHong Zhang     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
33fcd7ac73SHong Zhang   } else {
34fcd7ac73SHong Zhang     PetscBool assembled;
35fcd7ac73SHong Zhang     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
36fcd7ac73SHong Zhang     if (assembled) {
37fcd7ac73SHong Zhang       ierr = MatZeroEntries(J);CHKERRQ(ierr);
38fcd7ac73SHong Zhang     }
39fcd7ac73SHong Zhang   }
40fcd7ac73SHong Zhang 
41*f6af9589SHong Zhang   if (!coloring->vscale) { //ctype == IS_COLORING_GHOSTED ?
42*f6af9589SHong Zhang     printf("[%d] create a non-ghostedcoloring->vscale\n",rank);
43*f6af9589SHong Zhang     ierr = VecDuplicate(x1,&coloring->vscale);CHKERRQ(ierr);
44*f6af9589SHong Zhang     //SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->scale");
45fcd7ac73SHong Zhang   }
46fcd7ac73SHong Zhang 
47fcd7ac73SHong Zhang   /* Set w1 = F(x1) */
48fcd7ac73SHong Zhang   if (!coloring->fset) {
49fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
50*f6af9589SHong Zhang     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
51fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
52fcd7ac73SHong Zhang   } else {
53fcd7ac73SHong Zhang     coloring->fset = PETSC_FALSE;
54fcd7ac73SHong Zhang   }
55fcd7ac73SHong Zhang 
56fcd7ac73SHong Zhang   if (!coloring->w3) {
57*f6af9589SHong Zhang     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
58fcd7ac73SHong Zhang     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
59fcd7ac73SHong Zhang   }
60fcd7ac73SHong Zhang   w3 = coloring->w3;
61fcd7ac73SHong Zhang 
62*f6af9589SHong Zhang   /* Compute vscale = dx - the local scale factors, including ghost points */
63*f6af9589SHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* =J's row ownershipRange, used by ghosted x! */
64*f6af9589SHong Zhang   //printf("[%d] start end %d %d\n",rank,start,end);
65*f6af9589SHong Zhang 
66*f6af9589SHong Zhang   if (coloring->htype[0] == 'w') {
67*f6af9589SHong Zhang     /* tacky test; need to make systematic if we add other approaches to computing h*/
68*f6af9589SHong Zhang     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
69*f6af9589SHong Zhang     dx = (1.0 + unorm)*epsilon;
70*f6af9589SHong Zhang     ierr = VecSet(coloring->vscale,dx);CHKERRQ(ierr);
71*f6af9589SHong Zhang   } else { // buggy!!!
72*f6af9589SHong Zhang     ierr = VecGetLocalSize(x1,&N);CHKERRQ(ierr);
73*f6af9589SHong Zhang     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
74fcd7ac73SHong Zhang     ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
75fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GHOSTED) {
76fcd7ac73SHong Zhang       col_start = 0; col_end = N;
77fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GLOBAL) {
78*f6af9589SHong Zhang       xx           = xx - start;  //???
79*f6af9589SHong Zhang       vscale_array = vscale_array - start; //???
80fcd7ac73SHong Zhang       col_start    = start; col_end = N + start;
81fcd7ac73SHong Zhang     }
82fcd7ac73SHong Zhang     for (col=col_start; col<col_end; col++) {
83fcd7ac73SHong Zhang       /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
84fcd7ac73SHong Zhang       dx = xx[col];
85fcd7ac73SHong Zhang       if (dx == (PetscScalar)0.0) dx = 1.0;
86fcd7ac73SHong Zhang       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
87fcd7ac73SHong Zhang       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
88fcd7ac73SHong Zhang       dx               *= epsilon;
89*f6af9589SHong Zhang       vscale_array[col] = (PetscScalar)dx;
90*f6af9589SHong Zhang     }
91fcd7ac73SHong Zhang   }
92fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) vscale_array = vscale_array + start;
93fcd7ac73SHong Zhang   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
94fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) {
95fcd7ac73SHong Zhang     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
96fcd7ac73SHong Zhang     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
97fcd7ac73SHong Zhang   }
98fcd7ac73SHong Zhang 
99fcd7ac73SHong Zhang   if (coloring->vscaleforrow) {
100fcd7ac73SHong Zhang     vscaleforrow = coloring->vscaleforrow;
101fcd7ac73SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
102fcd7ac73SHong Zhang 
1031b97d346SHong Zhang   /* Loop over each color */
104f5aae955SHong Zhang   nz = 0;
105fcd7ac73SHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
106fcd7ac73SHong Zhang   for (k=0; k<coloring->ncolors; k++) {
107fcd7ac73SHong Zhang     coloring->currentcolor = k;
108fcd7ac73SHong Zhang 
109fcd7ac73SHong Zhang     /*
110fcd7ac73SHong Zhang       Loop over each column associated with color
111*f6af9589SHong Zhang       adding the perturbation to the vector w3 = x1 + dx.
112fcd7ac73SHong Zhang     */
113*f6af9589SHong Zhang     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
114*f6af9589SHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
115*f6af9589SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
116fcd7ac73SHong Zhang     for (l=0; l<coloring->ncolumns[k]; l++) {
117*f6af9589SHong Zhang       col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
118*f6af9589SHong Zhang       dx             = vscale_array[coloring->vscaleforrow[k][l]];
119fcd7ac73SHong Zhang       w3_array[col] += dx;
120*f6af9589SHong Zhang       /*
121*f6af9589SHong Zhang       if (rank ==1) {
122*f6af9589SHong Zhang         printf("[%d] w3[%d] += vscale[%d] %g\n",rank,col,coloring->vscaleforrow[k][l],dx);
123*f6af9589SHong Zhang       }
124*f6af9589SHong Zhang        */
125fcd7ac73SHong Zhang     }
126fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
127fcd7ac73SHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
128fcd7ac73SHong Zhang 
129fcd7ac73SHong Zhang     /*
130fcd7ac73SHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
131fcd7ac73SHong Zhang                            w2 = F(x1 + dx) - F(x1)
132fcd7ac73SHong Zhang     */
133fcd7ac73SHong Zhang     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
134fcd7ac73SHong Zhang     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
135fcd7ac73SHong Zhang     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
136fcd7ac73SHong Zhang     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
137fcd7ac73SHong Zhang 
1381b97d346SHong Zhang     /* Loop over rows of vector, putting results into Jacobian matrix */
139fcd7ac73SHong Zhang     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
140fcd7ac73SHong Zhang     for (l=0; l<coloring->nrows[k]; l++) {
141fcd7ac73SHong Zhang       row                        = coloring->rows[k][l];   /* local row index */
142*f6af9589SHong Zhang       *(coloring->valaddr[nz++]) = y[row]/vscale_array[vscaleforrow[k][l]];
143fcd7ac73SHong Zhang     }
144fcd7ac73SHong Zhang     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
145fcd7ac73SHong Zhang   } /* endof for each color */
146f5aae955SHong 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);
147f5aae955SHong Zhang   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148f5aae955SHong Zhang   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149fcd7ac73SHong Zhang 
150fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
151fcd7ac73SHong Zhang   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
152*f6af9589SHong Zhang   ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
153fcd7ac73SHong Zhang 
154fcd7ac73SHong Zhang   coloring->currentcolor = -1;
155fcd7ac73SHong Zhang   PetscFunctionReturn(0);
156fcd7ac73SHong Zhang }
157fcd7ac73SHong Zhang 
158fcd7ac73SHong Zhang #undef __FUNCT__
159fcd7ac73SHong Zhang #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new"
160fcd7ac73SHong Zhang PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c)
161fcd7ac73SHong Zhang {
162fcd7ac73SHong Zhang   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
163fcd7ac73SHong Zhang   PetscErrorCode         ierr;
164fcd7ac73SHong Zhang   PetscMPIInt            size,*ncolsonproc,*disp,nn;
165fcd7ac73SHong Zhang   PetscInt               i,n,nrows,j,k,m,ncols,col;
166d3825b63SHong Zhang   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog;
167fcd7ac73SHong Zhang   PetscInt               nis = iscoloring->n,nctot,*cols;
168d3825b63SHong Zhang   PetscInt               *rowhit,cstart,cend,colb;
169fcd7ac73SHong Zhang   PetscInt               *columnsforrow,l;
170fcd7ac73SHong Zhang   IS                     *isa;
171fcd7ac73SHong Zhang   ISLocalToGlobalMapping map = mat->cmap->mapping;
172fcd7ac73SHong Zhang   PetscInt               ctype=c->ctype;
173fcd7ac73SHong Zhang   Mat                    A=aij->A,B=aij->B;
174fcd7ac73SHong Zhang   Mat_SeqAIJ             *cspA=(Mat_SeqAIJ*)A->data, *cspB=(Mat_SeqAIJ*)B->data;
175fcd7ac73SHong Zhang   PetscScalar            *A_val=cspA->a,*B_val=cspB->a;
176fcd7ac73SHong Zhang   PetscInt               spidx;
1771b97d346SHong Zhang   PetscInt               *spidxA,*spidxB,nz;
178d3825b63SHong Zhang   PetscScalar            **valaddrhit;
179fcd7ac73SHong Zhang 
180fcd7ac73SHong Zhang   PetscFunctionBegin;
181fcd7ac73SHong 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");
182fcd7ac73SHong Zhang 
183fcd7ac73SHong Zhang   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
184fcd7ac73SHong Zhang   else     ltog = NULL;
185fcd7ac73SHong Zhang   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
186fcd7ac73SHong Zhang 
187d3825b63SHong Zhang   m         = mat->rmap->n;
188fcd7ac73SHong Zhang   cstart    = mat->cmap->rstart;
189fcd7ac73SHong Zhang   cend      = mat->cmap->rend;
190fcd7ac73SHong Zhang   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
191fcd7ac73SHong Zhang   c->N      = mat->cmap->N;
192d3825b63SHong Zhang   c->m      = m;
193fcd7ac73SHong Zhang   c->rstart = mat->rmap->rstart;
194fcd7ac73SHong Zhang 
195fcd7ac73SHong Zhang   c->ncolors = nis;
196fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
197fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
198fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
199fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
200fcd7ac73SHong Zhang   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
201fcd7ac73SHong Zhang   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
202fcd7ac73SHong Zhang 
203d3825b63SHong Zhang   /* Allow access to data structures of local part of matrix
204d3825b63SHong Zhang    - creates aij->colmap which maps global column number to local number in part B */
205fcd7ac73SHong Zhang   if (!aij->colmap) {
206fcd7ac73SHong Zhang     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
207fcd7ac73SHong Zhang   }
208fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
209fcd7ac73SHong Zhang   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
210fcd7ac73SHong Zhang 
211d3825b63SHong Zhang   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
212d3825b63SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
213fcd7ac73SHong Zhang   nz = cspA->nz + cspB->nz; /* total nonzero entries of mat */
214d3825b63SHong Zhang   ierr = PetscMalloc(nz*sizeof(PetscScalar*),&c->valaddr);CHKERRQ(ierr);
215fcd7ac73SHong Zhang 
216fcd7ac73SHong Zhang   nz = 0;
217d3825b63SHong Zhang   for (i=0; i<nis; i++) { /* for each local color */
218fcd7ac73SHong Zhang     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
219fcd7ac73SHong Zhang     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
220fcd7ac73SHong Zhang 
221fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
222fcd7ac73SHong Zhang     if (n) {
223fcd7ac73SHong Zhang       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
224fcd7ac73SHong Zhang       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
225fcd7ac73SHong Zhang       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
226fcd7ac73SHong Zhang     } else {
227fcd7ac73SHong Zhang       c->columns[i] = 0;
228fcd7ac73SHong Zhang     }
229fcd7ac73SHong Zhang 
230fcd7ac73SHong Zhang     if (ctype == IS_COLORING_GLOBAL) {
231d3825b63SHong Zhang       /* Determine nctot, the total (parallel) number of columns of this color */
232fcd7ac73SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
233fcd7ac73SHong Zhang       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
234fcd7ac73SHong Zhang 
235d3825b63SHong Zhang       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
236fcd7ac73SHong Zhang       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
237fcd7ac73SHong Zhang       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
238fcd7ac73SHong Zhang       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
239fcd7ac73SHong Zhang       if (!nctot) {
240fcd7ac73SHong Zhang         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
241fcd7ac73SHong Zhang       }
242fcd7ac73SHong Zhang 
243fcd7ac73SHong Zhang       disp[0] = 0;
244fcd7ac73SHong Zhang       for (j=1; j<size; j++) {
245fcd7ac73SHong Zhang         disp[j] = disp[j-1] + ncolsonproc[j-1];
246fcd7ac73SHong Zhang       }
247fcd7ac73SHong Zhang 
248d3825b63SHong Zhang       /* Get cols, the complete list of columns for this color on each process */
249fcd7ac73SHong Zhang       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
250fcd7ac73SHong Zhang       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
251fcd7ac73SHong Zhang       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
252fcd7ac73SHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
253fcd7ac73SHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
254fcd7ac73SHong Zhang       nctot = n;
255fcd7ac73SHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
256fcd7ac73SHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
257fcd7ac73SHong Zhang     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
258fcd7ac73SHong Zhang 
2591b97d346SHong Zhang     /* Mark all rows affect by these columns */
260d3825b63SHong Zhang     ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
2611b97d346SHong Zhang     for (j=0; j<nctot; j++) { /* loop over columns*/
262fcd7ac73SHong Zhang       if (ctype == IS_COLORING_GHOSTED) {
263fcd7ac73SHong Zhang         col = ltog[cols[j]];
264fcd7ac73SHong Zhang       } else {
265fcd7ac73SHong Zhang         col = cols[j];
266fcd7ac73SHong Zhang       }
267fcd7ac73SHong Zhang       if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */
268d3825b63SHong Zhang         row = A_cj + A_ci[col-cstart];
269d3825b63SHong Zhang         nrows = A_ci[col-cstart+1] - A_ci[col-cstart];
270fcd7ac73SHong Zhang         /* loop over columns of A marking them in rowhit */
271d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
272d3825b63SHong Zhang           /* set valaddrhit for part A */
273fcd7ac73SHong Zhang           spidx            = spidxA[A_ci[col-cstart] + k];
274d3825b63SHong Zhang           valaddrhit[*row] = &A_val[spidx];
275d3825b63SHong Zhang           rowhit[*row++]   = col + 1;
276fcd7ac73SHong Zhang         }
277fcd7ac73SHong Zhang       } else { /* column is in off-diagonal block of matrix B */
278fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
279fcd7ac73SHong Zhang         ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
280fcd7ac73SHong Zhang         colb--;
281fcd7ac73SHong Zhang #else
282fcd7ac73SHong Zhang         colb = aij->colmap[col] - 1; /* local column index */
283fcd7ac73SHong Zhang #endif
284fcd7ac73SHong Zhang         if (colb == -1) {
285d3825b63SHong Zhang           nrows = 0;
286fcd7ac73SHong Zhang         } else {
287d3825b63SHong Zhang           row   = B_cj + B_ci[colb];
288d3825b63SHong Zhang           nrows = B_ci[colb+1] - B_ci[colb];
289fcd7ac73SHong Zhang         }
290fcd7ac73SHong Zhang         /* loop over columns of B marking them in rowhit */
291d3825b63SHong Zhang         for (k=0; k<nrows; k++) {
292d3825b63SHong Zhang           /* set valaddrhit for part B */
293fcd7ac73SHong Zhang           spidx            = spidxB[B_ci[colb] + k];
294d3825b63SHong Zhang           valaddrhit[*row] = &B_val[spidx];
295d3825b63SHong Zhang           rowhit[*row++] = col + 1;
296fcd7ac73SHong Zhang         }
297fcd7ac73SHong Zhang       }
298fcd7ac73SHong Zhang     }
299fcd7ac73SHong Zhang 
300fcd7ac73SHong Zhang     /* count the number of hits */
301fcd7ac73SHong Zhang     nrows = 0;
302d3825b63SHong Zhang     for (j=0; j<m; j++) {
303fcd7ac73SHong Zhang       if (rowhit[j]) nrows++;
304fcd7ac73SHong Zhang     }
305fcd7ac73SHong Zhang     c->nrows[i] = nrows;
306fcd7ac73SHong Zhang     ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
307fcd7ac73SHong Zhang     ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
308fcd7ac73SHong Zhang     ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
309fcd7ac73SHong Zhang     nrows       = 0;
310d3825b63SHong Zhang     for (j=0; j<m; j++) {
311fcd7ac73SHong Zhang       if (rowhit[j]) {
312fcd7ac73SHong Zhang         c->rows[i][nrows]          = j;              /* local row index */
313fcd7ac73SHong Zhang         c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
314d3825b63SHong Zhang         c->valaddr[nz++]           = valaddrhit[j];   /* address of mat value for this entry */
315fcd7ac73SHong Zhang         nrows++;
316fcd7ac73SHong Zhang       }
317fcd7ac73SHong Zhang     }
318fcd7ac73SHong Zhang     ierr = PetscFree(cols);CHKERRQ(ierr);
319fcd7ac73SHong Zhang   }
320fcd7ac73SHong 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);
321fcd7ac73SHong Zhang 
322d3825b63SHong Zhang   /* Optimize by adding the vscale, and scaleforrow[][] fields;
323*f6af9589SHong Zhang      vscale contains the "diagonal" on processor scalings followed by the off processor */
324fcd7ac73SHong Zhang   if (ctype == IS_COLORING_GLOBAL) {
325d3825b63SHong Zhang     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),m,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
326fcd7ac73SHong Zhang     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
327fcd7ac73SHong Zhang     for (k=0; k<c->ncolors; k++) {
328fcd7ac73SHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
329fcd7ac73SHong Zhang       for (l=0; l<c->nrows[k]; l++) {
330fcd7ac73SHong Zhang         col = c->columnsforrow[k][l];
331fcd7ac73SHong Zhang         if (col >= cstart && col < cend) {
332fcd7ac73SHong Zhang           /* column is in diagonal block of matrix */
333fcd7ac73SHong Zhang           colb = col - cstart;
334fcd7ac73SHong Zhang         } else {
335fcd7ac73SHong Zhang           /* column  is in "off-processor" part */
336fcd7ac73SHong Zhang #if defined(PETSC_USE_CTABLE)
337fcd7ac73SHong Zhang           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
338fcd7ac73SHong Zhang           colb--;
339fcd7ac73SHong Zhang #else
340fcd7ac73SHong Zhang           colb = aij->colmap[col] - 1;
341fcd7ac73SHong Zhang #endif
342fcd7ac73SHong Zhang           colb += cend - cstart;
343fcd7ac73SHong Zhang         }
344d3825b63SHong Zhang         c->vscaleforrow[k][l] = colb; /* local column index */
345fcd7ac73SHong Zhang       }
346fcd7ac73SHong Zhang     }
347fcd7ac73SHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
348fcd7ac73SHong Zhang     /* Get gtol mapping */
349fcd7ac73SHong Zhang     PetscInt N = mat->cmap->N,nlocal,*gtol;
350fcd7ac73SHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
351fcd7ac73SHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
352fcd7ac73SHong Zhang     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
353fcd7ac73SHong Zhang     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
354fcd7ac73SHong Zhang 
355fcd7ac73SHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
356fcd7ac73SHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
357fcd7ac73SHong Zhang     for (k=0; k<c->ncolors; k++) {
358fcd7ac73SHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
359fcd7ac73SHong Zhang       for (l=0; l<c->nrows[k]; l++) {
360fcd7ac73SHong Zhang         col                   = c->columnsforrow[k][l]; /* global column index */
361fcd7ac73SHong Zhang         c->vscaleforrow[k][l] = gtol[col];              /* local column index */
362fcd7ac73SHong Zhang       }
363fcd7ac73SHong Zhang     }
364fcd7ac73SHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
365fcd7ac73SHong Zhang   }
366fcd7ac73SHong Zhang   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
367fcd7ac73SHong Zhang 
368d3825b63SHong Zhang   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
369fcd7ac73SHong Zhang   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
370fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
371fcd7ac73SHong Zhang   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
372fcd7ac73SHong Zhang   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
373fcd7ac73SHong Zhang 
374fcd7ac73SHong Zhang   mat->ops->fdcoloringapply = MatFDColoringApply_MPIAIJ;
375fcd7ac73SHong Zhang   PetscFunctionReturn(0);
376fcd7ac73SHong Zhang }
377fcd7ac73SHong Zhang 
378fcd7ac73SHong Zhang /*------------------------------------------------------*/
379fcd7ac73SHong Zhang #undef __FUNCT__
3804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
381dfbe8321SBarry Smith PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
382a64fbb32SBarry Smith {
3836eaac0f3SBarry Smith   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
3846849ba73SBarry Smith   PetscErrorCode         ierr;
385b1d57f15SBarry Smith   PetscMPIInt            size,*ncolsonproc,*disp,nn;
3861a83f524SJed Brown   PetscInt               i,n,nrows,j,k,m,ncols,col;
387afcb2eb5SJed Brown   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
3881a83f524SJed Brown   PetscInt               nis = iscoloring->n,nctot,*cols;
389f6d58c54SBarry Smith   PetscInt               *rowhit,M,cstart,cend,colb;
390b1d57f15SBarry Smith   PetscInt               *columnsforrow,l;
391b9617806SBarry Smith   IS                     *isa;
392ace3abfcSBarry Smith   PetscBool              done,flg;
393992144d0SBarry Smith   ISLocalToGlobalMapping map   = mat->cmap->mapping;
394afcb2eb5SJed Brown   PetscInt               ctype=c->ctype;
395fcd7ac73SHong Zhang   PetscBool              new_impl=PETSC_FALSE;
396a64fbb32SBarry Smith 
3973a40ed3dSBarry Smith   PetscFunctionBegin;
398fcd7ac73SHong Zhang   ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr);
399fcd7ac73SHong Zhang   if (new_impl){
400fcd7ac73SHong Zhang     ierr =  MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr);
401fcd7ac73SHong Zhang     PetscFunctionReturn(0);
402fcd7ac73SHong Zhang   }
403ce94432eSBarry 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");
404522c5e43SBarry Smith 
405afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
406afcb2eb5SJed Brown   else     ltog = NULL;
407b9617806SBarry Smith   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
4083acb8795SBarry Smith 
409f6d58c54SBarry Smith   M         = mat->rmap->n;
410f6d58c54SBarry Smith   cstart    = mat->cmap->rstart;
411f6d58c54SBarry Smith   cend      = mat->cmap->rend;
412f6d58c54SBarry Smith   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
413f6d58c54SBarry Smith   c->N      = mat->cmap->N;
414f6d58c54SBarry Smith   c->m      = mat->rmap->n;
415f6d58c54SBarry Smith   c->rstart = mat->rmap->rstart;
416005c665bSBarry Smith 
417a64fbb32SBarry Smith   c->ncolors = nis;
418b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
419b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
420b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
421b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
422b1d57f15SBarry Smith   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
4233bb1ff40SBarry Smith   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
4246eaac0f3SBarry Smith 
4256eaac0f3SBarry Smith   /* Allow access to data structures of local part of matrix */
4266eaac0f3SBarry Smith   if (!aij->colmap) {
427ab9863d7SBarry Smith     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
4286eaac0f3SBarry Smith   }
4293acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
4303acb8795SBarry Smith   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
4316eaac0f3SBarry Smith 
432b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
433b1d57f15SBarry Smith   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
4346eaac0f3SBarry Smith 
435a64fbb32SBarry Smith   for (i=0; i<nis; i++) {
436b9b97703SBarry Smith     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
437a64fbb32SBarry Smith     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
4382205254eSKarl Rupp 
439fcd7ac73SHong Zhang     c->ncolumns[i] = n; /* local number of columns of this color on this process */
440a64fbb32SBarry Smith     if (n) {
441b1d57f15SBarry Smith       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
4423bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
443b1d57f15SBarry Smith       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
444a64fbb32SBarry Smith     } else {
445a64fbb32SBarry Smith       c->columns[i] = 0;
446a64fbb32SBarry Smith     }
447a64fbb32SBarry Smith 
4488ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
4496eaac0f3SBarry Smith       /* Determine the total (parallel) number of columns of this color */
450ce94432eSBarry Smith       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
451687f1162SBarry Smith       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
452b8f8c88eSHong Zhang 
4534dc2109aSBarry Smith       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
454ce94432eSBarry Smith       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4552205254eSKarl Rupp       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
4563a7fca6bSBarry Smith       if (!nctot) {
457ae15b995SBarry Smith         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
4583a7fca6bSBarry Smith       }
4596eaac0f3SBarry Smith 
4606eaac0f3SBarry Smith       disp[0] = 0;
4616eaac0f3SBarry Smith       for (j=1; j<size; j++) {
4626eaac0f3SBarry Smith         disp[j] = disp[j-1] + ncolsonproc[j-1];
4636eaac0f3SBarry Smith       }
4646eaac0f3SBarry Smith 
4656eaac0f3SBarry Smith       /* Get complete list of columns for color on each processor */
466b1d57f15SBarry Smith       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
467ce94432eSBarry Smith       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
4681d79065fSBarry Smith       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
469b8f8c88eSHong Zhang     } else if (ctype == IS_COLORING_GHOSTED) {
470b8f8c88eSHong Zhang       /* Determine local number of columns of this color on this process, including ghost points */
471b8f8c88eSHong Zhang       nctot = n;
472b8f8c88eSHong Zhang       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
473b8f8c88eSHong Zhang       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
474f23aa3ddSBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
4756eaac0f3SBarry Smith 
4766eaac0f3SBarry Smith     /*
4776eaac0f3SBarry Smith        Mark all rows affect by these columns
4786eaac0f3SBarry Smith     */
479b8f8c88eSHong Zhang     /* Temporary option to allow for debugging/testing */
48090d69ab7SBarry Smith     flg  = PETSC_FALSE;
4810298fd71SBarry Smith     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
482f158e583SBarry Smith     if (!flg) { /*-----------------------------------------------------------------------------*/
483f158e583SBarry Smith       /* crude, fast version */
484b1d57f15SBarry Smith       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
485a64fbb32SBarry Smith       /* loop over columns*/
4866eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
487b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
488b8f8c88eSHong Zhang           col = ltog[cols[j]];
489b8f8c88eSHong Zhang         } else {
4906eaac0f3SBarry Smith           col = cols[j];
491b8f8c88eSHong Zhang         }
4926eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
4936eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
4946eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
4956eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
4966eaac0f3SBarry Smith         } else {
497aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
498cb9801acSJed Brown           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
499fa46199cSSatish Balay           colb--;
500b3d2dc96SSatish Balay #else
5016eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
502b3d2dc96SSatish Balay #endif
5036eaac0f3SBarry Smith           if (colb == -1) {
5046eaac0f3SBarry Smith             m = 0;
5056eaac0f3SBarry Smith           } else {
5066eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
5076eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
5086eaac0f3SBarry Smith           }
5096eaac0f3SBarry Smith         }
510a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
511a64fbb32SBarry Smith         for (k=0; k<m; k++) {
512a64fbb32SBarry Smith           rowhit[*rows++] = col + 1;
513a64fbb32SBarry Smith         }
514a64fbb32SBarry Smith       }
5156eaac0f3SBarry Smith 
516a64fbb32SBarry Smith       /* count the number of hits */
517a64fbb32SBarry Smith       nrows = 0;
5186eaac0f3SBarry Smith       for (j=0; j<M; j++) {
519a64fbb32SBarry Smith         if (rowhit[j]) nrows++;
520a64fbb32SBarry Smith       }
521a64fbb32SBarry Smith       c->nrows[i] = nrows;
522b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
523b1d57f15SBarry Smith       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
5243bb1ff40SBarry Smith       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
525a64fbb32SBarry Smith       nrows       = 0;
5266eaac0f3SBarry Smith       for (j=0; j<M; j++) {
527a64fbb32SBarry Smith         if (rowhit[j]) {
528fcd7ac73SHong Zhang           c->rows[i][nrows]          = j;              /* local row index */
529fcd7ac73SHong Zhang           c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
530a64fbb32SBarry Smith           nrows++;
531a64fbb32SBarry Smith         }
532a64fbb32SBarry Smith       }
533a64fbb32SBarry Smith     } else { /*-------------------------------------------------------------------------------*/
534f158e583SBarry Smith       /* slow version, using rowhit as a linked list */
535b1d57f15SBarry Smith       PetscInt currentcol,fm,mfm;
5366eaac0f3SBarry Smith       rowhit[M] = M;
537a64fbb32SBarry Smith       nrows     = 0;
538a64fbb32SBarry Smith       /* loop over columns*/
5396eaac0f3SBarry Smith       for (j=0; j<nctot; j++) {
540b8f8c88eSHong Zhang         if (ctype == IS_COLORING_GHOSTED) {
541b8f8c88eSHong Zhang           col = ltog[cols[j]];
542b8f8c88eSHong Zhang         } else {
5436eaac0f3SBarry Smith           col = cols[j];
544b8f8c88eSHong Zhang         }
5456eaac0f3SBarry Smith         if (col >= cstart && col < cend) {
5466eaac0f3SBarry Smith           /* column is in diagonal block of matrix */
5476eaac0f3SBarry Smith           rows = A_cj + A_ci[col-cstart];
5486eaac0f3SBarry Smith           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
5496eaac0f3SBarry Smith         } else {
550aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
5510f5bd95cSBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
552fa46199cSSatish Balay           colb--;
553b3d2dc96SSatish Balay #else
5546eaac0f3SBarry Smith           colb = aij->colmap[col] - 1;
555b3d2dc96SSatish Balay #endif
5566eaac0f3SBarry Smith           if (colb == -1) {
5576eaac0f3SBarry Smith             m = 0;
5586eaac0f3SBarry Smith           } else {
5596eaac0f3SBarry Smith             rows = B_cj + B_ci[colb];
5606eaac0f3SBarry Smith             m    = B_ci[colb+1] - B_ci[colb];
5616eaac0f3SBarry Smith           }
5626eaac0f3SBarry Smith         }
563b8f8c88eSHong Zhang 
564a64fbb32SBarry Smith         /* loop over columns marking them in rowhit */
5656eaac0f3SBarry Smith         fm = M;    /* fm points to first entry in linked list */
566a64fbb32SBarry Smith         for (k=0; k<m; k++) {
567a64fbb32SBarry Smith           currentcol = *rows++;
568a64fbb32SBarry Smith           /* is it already in the list? */
569a64fbb32SBarry Smith           do {
570a64fbb32SBarry Smith             mfm = fm;
571a64fbb32SBarry Smith             fm  = rowhit[fm];
572a64fbb32SBarry Smith           } while (fm < currentcol);
573a64fbb32SBarry Smith           /* not in list so add it */
574a64fbb32SBarry Smith           if (fm != currentcol) {
575a64fbb32SBarry Smith             nrows++;
576a64fbb32SBarry Smith             columnsforrow[currentcol] = col;
577a64fbb32SBarry Smith             /* next three lines insert new entry into linked list */
578a64fbb32SBarry Smith             rowhit[mfm]        = currentcol;
579a64fbb32SBarry Smith             rowhit[currentcol] = fm;
580a64fbb32SBarry Smith             fm                 = currentcol;
581a64fbb32SBarry Smith             /* fm points to present position in list since we know the columns are sorted */
582f23aa3ddSBarry Smith           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
583a64fbb32SBarry Smith         }
584a64fbb32SBarry Smith       }
585a64fbb32SBarry Smith       c->nrows[i] = nrows;
5862205254eSKarl Rupp 
587b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
588b1d57f15SBarry Smith       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
5893bb1ff40SBarry Smith       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
590a64fbb32SBarry Smith       /* now store the linked list of rows into c->rows[i] */
591a64fbb32SBarry Smith       nrows = 0;
5926eaac0f3SBarry Smith       fm    = rowhit[M];
593a64fbb32SBarry Smith       do {
594a64fbb32SBarry Smith         c->rows[i][nrows]            = fm;
595a64fbb32SBarry Smith         c->columnsforrow[i][nrows++] = columnsforrow[fm];
596a64fbb32SBarry Smith         fm                           = rowhit[fm];
5976eaac0f3SBarry Smith       } while (fm < M);
5986eaac0f3SBarry Smith     } /* ---------------------------------------------------------------------------------------*/
599606d414cSSatish Balay     ierr = PetscFree(cols);CHKERRQ(ierr);
6006eaac0f3SBarry Smith   }
60130b34957SBarry Smith 
60230b34957SBarry Smith   /* Optimize by adding the vscale, and scaleforrow[][] fields */
60330b34957SBarry Smith   /*
60430b34957SBarry Smith        vscale will contain the "diagonal" on processor scalings followed by the off processor
60530b34957SBarry Smith   */
6068ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
607ce94432eSBarry Smith     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
608b1d57f15SBarry Smith     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
60930b34957SBarry Smith     for (k=0; k<c->ncolors; k++) {
610b1d57f15SBarry Smith       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
61130b34957SBarry Smith       for (l=0; l<c->nrows[k]; l++) {
61230b34957SBarry Smith         col = c->columnsforrow[k][l];
61330b34957SBarry Smith         if (col >= cstart && col < cend) {
61430b34957SBarry Smith           /* column is in diagonal block of matrix */
61530b34957SBarry Smith           colb = col - cstart;
61630b34957SBarry Smith         } else {
61730b34957SBarry Smith           /* column  is in "off-processor" part */
61830b34957SBarry Smith #if defined(PETSC_USE_CTABLE)
61930b34957SBarry Smith           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
62030b34957SBarry Smith           colb--;
62130b34957SBarry Smith #else
62230b34957SBarry Smith           colb = aij->colmap[col] - 1;
62330b34957SBarry Smith #endif
62430b34957SBarry Smith           colb += cend - cstart;
62530b34957SBarry Smith         }
62630b34957SBarry Smith         c->vscaleforrow[k][l] = colb;
62730b34957SBarry Smith       }
62830b34957SBarry Smith     }
629b8f8c88eSHong Zhang   } else if (ctype == IS_COLORING_GHOSTED) {
630b8f8c88eSHong Zhang     /* Get gtol mapping */
631afcb2eb5SJed Brown     PetscInt N = mat->cmap->N,nlocal,*gtol;
632b8f8c88eSHong Zhang     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
633b8f8c88eSHong Zhang     for (i=0; i<N; i++) gtol[i] = -1;
634afcb2eb5SJed Brown     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
635afcb2eb5SJed Brown     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
636b8f8c88eSHong Zhang 
637b8f8c88eSHong Zhang     c->vscale = 0; /* will be created in MatFDColoringApply() */
638b8f8c88eSHong Zhang     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
639b8f8c88eSHong Zhang     for (k=0; k<c->ncolors; k++) {
640b8f8c88eSHong Zhang       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
641b8f8c88eSHong Zhang       for (l=0; l<c->nrows[k]; l++) {
642b8f8c88eSHong Zhang         col = c->columnsforrow[k][l];      /* global column index */
6432205254eSKarl Rupp 
644b8f8c88eSHong Zhang         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
645b8f8c88eSHong Zhang       }
646b8f8c88eSHong Zhang     }
647b8f8c88eSHong Zhang     ierr = PetscFree(gtol);CHKERRQ(ierr);
648b8f8c88eSHong Zhang   }
649b9617806SBarry Smith   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
65030b34957SBarry Smith 
651606d414cSSatish Balay   ierr = PetscFree(rowhit);CHKERRQ(ierr);
652606d414cSSatish Balay   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
6533acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
6543acb8795SBarry Smith   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
655afcb2eb5SJed Brown   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
6563a40ed3dSBarry Smith   PetscFunctionReturn(0);
657a64fbb32SBarry Smith }
658a64fbb32SBarry Smith 
659b9617806SBarry Smith 
660b9617806SBarry Smith 
661b9617806SBarry Smith 
662b9617806SBarry Smith 
663b9617806SBarry Smith 
664