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