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