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