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