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