xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision c53567a07016831a3aeb692e1ce5d8d772267d49)
1 
2 #include <../src/mat/impls/aij/mpi/mpiaij.h>
3 
4 extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "MatFDColoringApply_AIJ_new"
8 PetscErrorCode  MatFDColoringApply_AIJ_new(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
9 {
10   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
11   PetscErrorCode ierr;
12   PetscInt       k,cstart,cend,l,row,col,nz;
13   PetscScalar    dx=0.0,*y,*xx,*w3_array;
14   PetscScalar    *vscale_array;
15   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
16   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
17   void           *fctx=coloring->fctx;
18   PetscBool      flg=PETSC_FALSE;
19   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
20   Mat_MPIAIJ     *aij=(Mat_MPIAIJ*)J->data;
21   MatEntry       *Jentry=coloring->matentry;
22   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
23 
24   PetscFunctionBegin;
25   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
26   ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr);
27   if (flg) {
28     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
29   } else {
30     PetscBool assembled;
31     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
32     if (assembled) {
33       ierr = MatZeroEntries(J);CHKERRQ(ierr);
34     }
35   }
36 
37   /* create vscale for storing dx */
38   if (!vscale) {
39     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
40       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
41     } else if (ctype == IS_COLORING_GHOSTED) {
42       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
43     }
44     coloring->vscale = vscale;
45   }
46 
47   /* (1) Set w1 = F(x1) */
48   if (!coloring->fset) {
49     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
50     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
51     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
52   } else {
53     coloring->fset = PETSC_FALSE;
54   }
55 
56   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
57   if (coloring->htype[0] == 'w') {
58     /* vscale = dx is a constant scalar */
59     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
60     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
61   } else {
62     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
63     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
64     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
65     for (col=0; col<nxloc; col++) {
66       dx = xx[col];
67       if (PetscAbsScalar(dx) < umin) {
68         if (PetscRealPart(dx) >= 0.0)      dx = umin;
69         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
70       }
71       dx               *= epsilon;
72       vscale_array[col] = 1.0/dx;
73     }
74     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
75     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
76   }
77   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') {
78     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
79     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
80   }
81 
82   /* (3) Loop over each color */
83   if (!coloring->w3) {
84     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
85     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
86   }
87   w3 = coloring->w3;
88 
89   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
90   if (vscale) {
91     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
92   }
93   nz   = 0;
94   for (k=0; k<ncolors; k++) {
95     coloring->currentcolor = k;
96 
97     /*
98       (3-1) Loop over each column associated with color
99       adding the perturbation to the vector w3 = x1 + dx.
100     */
101     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
102     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
103     if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
104     if (coloring->htype[0] == 'w') {
105       for (l=0; l<ncolumns[k]; l++) {
106         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
107         w3_array[col] += 1.0/dx;
108       }
109     } else { /* htype == 'ds' */
110       vscale_array -= cstart; /* shift pointer so global index can be used */
111       for (l=0; l<ncolumns[k]; l++) {
112         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
113         w3_array[col] += 1.0/vscale_array[col];
114       }
115       vscale_array += cstart;
116     }
117     if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
118     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
119 
120     /*
121       (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
122                            w2 = F(x1 + dx) - F(x1)
123     */
124     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
125     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
126     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
127     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
128 
129     /*
130      (3-3) Loop over rows of vector, putting results into Jacobian matrix
131     */
132     nrows_k = nrows[k];
133     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
134     if (coloring->htype[0] == 'w') {
135       for (l=0; l<nrows_k; l++) {
136         row                     = Jentry[nz].row;   /* local row index */
137         *(Jentry[nz++].valaddr) = y[row]*dx;
138       }
139     } else { /* htype == 'ds' */
140       for (l=0; l<nrows_k; l++) {
141         row                   = Jentry[nz].row;   /* local row index */
142         *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
143         nz++;
144       }
145     }
146     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
147   }
148   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150   if (vscale) {
151     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
152   }
153 
154   coloring->currentcolor = -1;
155   PetscFunctionReturn(0);
156 }
157 
158 #undef __FUNCT__
159 #define __FUNCT__ "MatFDColoringCreate_MPIAIJ_new"
160 PetscErrorCode MatFDColoringCreate_MPIAIJ_new(Mat mat,ISColoring iscoloring,MatFDColoring c)
161 {
162   Mat_MPIAIJ             *aij=(Mat_MPIAIJ*)mat->data;
163   PetscErrorCode         ierr;
164   PetscMPIInt            size,*ncolsonproc,*disp,nn;
165   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col;
166   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL;
167   PetscInt               nis=iscoloring->n,nctot,*cols;
168   PetscInt               *rowhit,cstart,cend,colb;
169   IS                     *isa;
170   ISLocalToGlobalMapping map=mat->cmap->mapping;
171   PetscInt               ctype=c->ctype;
172   Mat                    A=aij->A,B=aij->B;
173   Mat_SeqAIJ             *spA=(Mat_SeqAIJ*)A->data,*spB=(Mat_SeqAIJ*)B->data;
174   PetscScalar            *A_val=spA->a,*B_val=spB->a;
175   PetscInt               spidx;
176   PetscInt               *spidxA,*spidxB,nz;
177   PetscScalar            **valaddrhit;
178   MatEntry               *Jentry;
179 
180   PetscFunctionBegin;
181   if (ctype == IS_COLORING_GHOSTED) {
182     if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
183     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
184   }
185 
186   m         = mat->rmap->n;
187   cstart    = mat->cmap->rstart;
188   cend      = mat->cmap->rend;
189   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
190   c->N      = mat->cmap->N;
191   c->m      = m;
192   c->rstart = mat->rmap->rstart;
193 
194   c->ncolors = nis;
195   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
196   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
197   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
198   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
199 
200   nz         = spA->nz + spB->nz; /* total nonzero entries of mat */
201   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
202   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
203   c->matentry = Jentry;
204 
205   /* Allow access to data structures of local part of matrix
206    - creates aij->colmap which maps global column number to local number in part B */
207   if (!aij->colmap) {
208     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
209   }
210   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
211   ierr = MatGetColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
212 
213   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr);
214   nz = 0;
215   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
216   for (i=0; i<nis; i++) { /* for each local color */
217     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
218     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
219 
220     c->ncolumns[i] = n; /* local number of columns of this color on this process */
221     if (n) {
222       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
223       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
224       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
225       /* convert global column indices to local ones! */
226 
227     } else {
228       c->columns[i] = 0;
229     }
230 
231     if (ctype == IS_COLORING_GLOBAL) {
232       /* Determine nctot, the total (parallel) number of columns of this color */
233       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
234       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
235 
236       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
237       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
238       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
239       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
240       if (!nctot) {
241         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
242       }
243 
244       disp[0] = 0;
245       for (j=1; j<size; j++) {
246         disp[j] = disp[j-1] + ncolsonproc[j-1];
247       }
248 
249       /* Get cols, the complete list of columns for this color on each process */
250       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
251       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
252       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
253     } else if (ctype == IS_COLORING_GHOSTED) {
254       /* Determine local number of columns of this color on this process, including ghost points */
255       nctot = n;
256       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
257       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
258     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
259 
260     /* Mark all rows affect by these columns */
261     ierr = PetscMemzero(rowhit,m*sizeof(PetscInt));CHKERRQ(ierr);
262     nrows_i = 0;
263     for (j=0; j<nctot; j++) { /* loop over columns*/
264       if (ctype == IS_COLORING_GHOSTED) {
265         col = ltog[cols[j]];
266       } else {
267         col = cols[j];
268       }
269       if (col >= cstart && col < cend) { /* column is in diagonal block of matrix A */
270         row      = A_cj + A_ci[col-cstart];
271         nrows    = A_ci[col-cstart+1] - A_ci[col-cstart];
272         nrows_i += nrows;
273         /* loop over columns of A marking them in rowhit */
274         for (k=0; k<nrows; k++) {
275           /* set valaddrhit for part A */
276           spidx            = spidxA[A_ci[col-cstart] + k];
277           valaddrhit[*row] = &A_val[spidx];
278           rowhit[*row++]   = col - cstart + 1; /* local column index */
279         }
280       } else { /* column is in off-diagonal block of matrix B */
281 #if defined(PETSC_USE_CTABLE)
282         ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
283         colb--;
284 #else
285         colb = aij->colmap[col] - 1; /* local column index */
286 #endif
287         if (colb == -1) {
288           nrows = 0;
289         } else {
290           row   = B_cj + B_ci[colb];
291           nrows = B_ci[colb+1] - B_ci[colb];
292         }
293         nrows_i += nrows;
294         /* loop over columns of B marking them in rowhit */
295         for (k=0; k<nrows; k++) {
296           /* set valaddrhit for part B */
297           spidx            = spidxB[B_ci[colb] + k];
298           valaddrhit[*row] = &B_val[spidx];
299           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
300         }
301       }
302     }
303     c->nrows[i] = nrows_i;
304 
305     for (j=0; j<m; j++) {
306       if (rowhit[j]) {
307         Jentry[nz].row     = j;              /* local row index */
308         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
309         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
310         nz++;
311       }
312     }
313     ierr = PetscFree(cols);CHKERRQ(ierr);
314   }
315   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
316   if (nz != spA->nz + spB->nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"nz %d != mat->nz %d\n",nz,spA->nz+spB->nz);
317 
318   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
319   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
320   ierr = MatRestoreColumnIJ_SeqAIJ_Color(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
321   if (ctype == IS_COLORING_GHOSTED) {
322     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
323   }
324 
325   mat->ops->fdcoloringapply = MatFDColoringApply_AIJ_new;
326   PetscFunctionReturn(0);
327 }
328 
329 /*------------------------------------------------------*/
330 #undef __FUNCT__
331 #define __FUNCT__ "MatFDColoringCreate_MPIAIJ"
332 PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
333 {
334   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)mat->data;
335   PetscErrorCode         ierr;
336   PetscMPIInt            size,*ncolsonproc,*disp,nn;
337   PetscInt               i,n,nrows,j,k,m,ncols,col;
338   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*rows = 0,*ltog;
339   PetscInt               nis = iscoloring->n,nctot,*cols;
340   PetscInt               *rowhit,M,cstart,cend,colb;
341   PetscInt               *columnsforrow,l;
342   IS                     *isa;
343   PetscBool              done,flg;
344   ISLocalToGlobalMapping map   = mat->cmap->mapping;
345   PetscInt               ctype=c->ctype;
346   PetscBool              new_impl=PETSC_FALSE;
347 
348   PetscFunctionBegin;
349   ierr = PetscOptionsName("-new","using new impls","",&new_impl);CHKERRQ(ierr);
350   if (new_impl){
351     ierr =  MatFDColoringCreate_MPIAIJ_new(mat,iscoloring,c);CHKERRQ(ierr);
352     PetscFunctionReturn(0);
353   }
354   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");
355 
356   if (map) {ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);}
357   else     ltog = NULL;
358   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
359 
360   M         = mat->rmap->n;
361   cstart    = mat->cmap->rstart;
362   cend      = mat->cmap->rend;
363   c->M      = mat->rmap->N;         /* set the global rows and columns and local rows */
364   c->N      = mat->cmap->N;
365   c->m      = mat->rmap->n;
366   c->rstart = mat->rmap->rstart;
367 
368   c->ncolors = nis;
369   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
370   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
371   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
372   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->rows);CHKERRQ(ierr);
373   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columnsforrow);CHKERRQ(ierr);
374   ierr       = PetscLogObjectMemory((PetscObject)c,5*nis*sizeof(PetscInt));CHKERRQ(ierr);
375 
376   /* Allow access to data structures of local part of matrix */
377   if (!aij->colmap) {
378     ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
379   }
380   ierr = MatGetColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
381   ierr = MatGetColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
382 
383   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
384   ierr = PetscMalloc((M+1)*sizeof(PetscInt),&columnsforrow);CHKERRQ(ierr);
385 
386   for (i=0; i<nis; i++) {
387     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
388     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
389 
390     c->ncolumns[i] = n; /* local number of columns of this color on this process */
391     if (n) {
392       ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr);
393       ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr);
394       ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr);
395     } else {
396       c->columns[i] = 0;
397     }
398 
399     if (ctype == IS_COLORING_GLOBAL) {
400       /* Determine the total (parallel) number of columns of this color */
401       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr);
402       ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr);
403 
404       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
405       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
406       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
407       if (!nctot) {
408         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
409       }
410 
411       disp[0] = 0;
412       for (j=1; j<size; j++) {
413         disp[j] = disp[j-1] + ncolsonproc[j-1];
414       }
415 
416       /* Get complete list of columns for color on each processor */
417       ierr = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
418       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
419       ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
420     } else if (ctype == IS_COLORING_GHOSTED) {
421       /* Determine local number of columns of this color on this process, including ghost points */
422       nctot = n;
423       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&cols);CHKERRQ(ierr);
424       ierr  = PetscMemcpy(cols,is,n*sizeof(PetscInt));CHKERRQ(ierr);
425     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
426 
427     /*
428        Mark all rows affect by these columns
429     */
430     /* Temporary option to allow for debugging/testing */
431     flg  = PETSC_FALSE;
432     ierr = PetscOptionsGetBool(NULL,"-matfdcoloring_slow",&flg,NULL);CHKERRQ(ierr);
433     if (!flg) { /*-----------------------------------------------------------------------------*/
434       /* crude, fast version */
435       ierr = PetscMemzero(rowhit,M*sizeof(PetscInt));CHKERRQ(ierr);
436       /* loop over columns*/
437       for (j=0; j<nctot; j++) {
438         if (ctype == IS_COLORING_GHOSTED) {
439           col = ltog[cols[j]];
440         } else {
441           col = cols[j];
442         }
443         if (col >= cstart && col < cend) {
444           /* column is in diagonal block of matrix */
445           rows = A_cj + A_ci[col-cstart];
446           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
447         } else {
448 #if defined(PETSC_USE_CTABLE)
449           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
450           colb--;
451 #else
452           colb = aij->colmap[col] - 1;
453 #endif
454           if (colb == -1) {
455             m = 0;
456           } else {
457             rows = B_cj + B_ci[colb];
458             m    = B_ci[colb+1] - B_ci[colb];
459           }
460         }
461         /* loop over columns marking them in rowhit */
462         for (k=0; k<m; k++) {
463           rowhit[*rows++] = col + 1;
464         }
465       }
466 
467       /* count the number of hits */
468       nrows = 0;
469       for (j=0; j<M; j++) {
470         if (rowhit[j]) nrows++;
471       }
472       c->nrows[i] = nrows;
473       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
474       ierr        = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
475       ierr        = PetscLogObjectMemory((PetscObject)c,2*(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
476       nrows       = 0;
477       for (j=0; j<M; j++) {
478         if (rowhit[j]) {
479           c->rows[i][nrows]          = j;              /* local row index */
480           c->columnsforrow[i][nrows] = rowhit[j] - 1;  /* global column index */
481           nrows++;
482         }
483       }
484     } else { /*-------------------------------------------------------------------------------*/
485       /* slow version, using rowhit as a linked list */
486       PetscInt currentcol,fm,mfm;
487       rowhit[M] = M;
488       nrows     = 0;
489       /* loop over columns*/
490       for (j=0; j<nctot; j++) {
491         if (ctype == IS_COLORING_GHOSTED) {
492           col = ltog[cols[j]];
493         } else {
494           col = cols[j];
495         }
496         if (col >= cstart && col < cend) {
497           /* column is in diagonal block of matrix */
498           rows = A_cj + A_ci[col-cstart];
499           m    = A_ci[col-cstart+1] - A_ci[col-cstart];
500         } else {
501 #if defined(PETSC_USE_CTABLE)
502           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
503           colb--;
504 #else
505           colb = aij->colmap[col] - 1;
506 #endif
507           if (colb == -1) {
508             m = 0;
509           } else {
510             rows = B_cj + B_ci[colb];
511             m    = B_ci[colb+1] - B_ci[colb];
512           }
513         }
514 
515         /* loop over columns marking them in rowhit */
516         fm = M;    /* fm points to first entry in linked list */
517         for (k=0; k<m; k++) {
518           currentcol = *rows++;
519           /* is it already in the list? */
520           do {
521             mfm = fm;
522             fm  = rowhit[fm];
523           } while (fm < currentcol);
524           /* not in list so add it */
525           if (fm != currentcol) {
526             nrows++;
527             columnsforrow[currentcol] = col;
528             /* next three lines insert new entry into linked list */
529             rowhit[mfm]        = currentcol;
530             rowhit[currentcol] = fm;
531             fm                 = currentcol;
532             /* fm points to present position in list since we know the columns are sorted */
533           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid coloring of matrix detected");
534         }
535       }
536       c->nrows[i] = nrows;
537 
538       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->rows[i]);CHKERRQ(ierr);
539       ierr = PetscMalloc((nrows+1)*sizeof(PetscInt),&c->columnsforrow[i]);CHKERRQ(ierr);
540       ierr = PetscLogObjectMemory((PetscObject)c,(nrows+1)*sizeof(PetscInt));CHKERRQ(ierr);
541       /* now store the linked list of rows into c->rows[i] */
542       nrows = 0;
543       fm    = rowhit[M];
544       do {
545         c->rows[i][nrows]            = fm;
546         c->columnsforrow[i][nrows++] = columnsforrow[fm];
547         fm                           = rowhit[fm];
548       } while (fm < M);
549     } /* ---------------------------------------------------------------------------------------*/
550     ierr = PetscFree(cols);CHKERRQ(ierr);
551   }
552 
553   /* Optimize by adding the vscale, and scaleforrow[][] fields */
554   /*
555        vscale will contain the "diagonal" on processor scalings followed by the off processor
556   */
557   if (ctype == IS_COLORING_GLOBAL) {
558     ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),aij->A->rmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
559     ierr = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
560     for (k=0; k<c->ncolors; k++) {
561       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
562       for (l=0; l<c->nrows[k]; l++) {
563         col = c->columnsforrow[k][l];
564         if (col >= cstart && col < cend) {
565           /* column is in diagonal block of matrix */
566           colb = col - cstart;
567         } else {
568           /* column  is in "off-processor" part */
569 #if defined(PETSC_USE_CTABLE)
570           ierr = PetscTableFind(aij->colmap,col+1,&colb);CHKERRQ(ierr);
571           colb--;
572 #else
573           colb = aij->colmap[col] - 1;
574 #endif
575           colb += cend - cstart;
576         }
577         c->vscaleforrow[k][l] = colb;
578       }
579     }
580   } else if (ctype == IS_COLORING_GHOSTED) {
581     /* Get gtol mapping */
582     PetscInt N = mat->cmap->N,nlocal,*gtol;
583     ierr = PetscMalloc((N+1)*sizeof(PetscInt),&gtol);CHKERRQ(ierr);
584     for (i=0; i<N; i++) gtol[i] = -1;
585     ierr = ISLocalToGlobalMappingGetSize(map,&nlocal);CHKERRQ(ierr);
586     for (i=0; i<nlocal; i++) gtol[ltog[i]] = i;
587 
588     c->vscale = 0; /* will be created in MatFDColoringApply() */
589     ierr      = PetscMalloc(c->ncolors*sizeof(PetscInt*),&c->vscaleforrow);CHKERRQ(ierr);
590     for (k=0; k<c->ncolors; k++) {
591       ierr = PetscMalloc((c->nrows[k]+1)*sizeof(PetscInt),&c->vscaleforrow[k]);CHKERRQ(ierr);
592       for (l=0; l<c->nrows[k]; l++) {
593         col = c->columnsforrow[k][l];      /* global column index */
594         c->vscaleforrow[k][l] = gtol[col]; /* local column index */
595       }
596     }
597     ierr = PetscFree(gtol);CHKERRQ(ierr);
598   }
599   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
600 
601   ierr = PetscFree(rowhit);CHKERRQ(ierr);
602   ierr = PetscFree(columnsforrow);CHKERRQ(ierr);
603   ierr = MatRestoreColumnIJ(aij->A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&done);CHKERRQ(ierr);
604   ierr = MatRestoreColumnIJ(aij->B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&done);CHKERRQ(ierr);
605   if (map) {ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);}
606   PetscFunctionReturn(0);
607 }
608 
609 
610 
611 
612 
613 
614