xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision d880da652b9b0febbb366ced8120a3620a856c60)
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   /* create vscale for storing dx */
25   if (!vscale) {
26     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
27       /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */
28       Mat_MPIBAIJ    *baij=(Mat_MPIBAIJ*)J->data;
29       PetscInt       *garray;
30       ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr);
31       for (i=0; i<baij->B->cmap->n/bs; i++) {
32         for (j=0; j<bs; j++) {
33           garray[i*bs+j] = bs*baij->garray[i]+j;
34         }
35       }
36       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr);
37       ierr = PetscFree(garray);CHKERRQ(ierr);
38     } else if (ctype == IS_COLORING_GHOSTED) {
39       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
40     }
41     coloring->vscale = vscale;
42   }
43 
44   /* (1) Set w1 = F(x1) */
45   if (!coloring->fset) {
46     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
47     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
48     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
49   } else {
50     coloring->fset = PETSC_FALSE;
51   }
52 
53   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
54   ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
55   if (coloring->htype[0] == 'w') {
56     /* vscale = dx is a constant scalar */
57     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
58     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
59   } else {
60     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
61     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
62     for (col=0; col<nxloc; col++) {
63       dx = xx[col];
64       if (PetscAbsScalar(dx) < umin) {
65         if (PetscRealPart(dx) >= 0.0)      dx = umin;
66         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
67       }
68       dx               *= epsilon;
69       vscale_array[col] = 1.0/dx;
70     }
71     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
72     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
73   }
74   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
75     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
76     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
77   }
78 
79   /* (3) Loop over each color */
80   if (!coloring->w3) {
81     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
82     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
83   }
84   w3 = coloring->w3;
85 
86   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
87   if (vscale) {
88     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
89   }
90   nz   = 0;
91   for (k=0; k<ncolors; k++) {
92     coloring->currentcolor = k;
93 
94     /*
95       (3-1) Loop over each column associated with color
96       adding the perturbation to the vector w3 = x1 + dx.
97     */
98     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
99     dy_i = dy;
100     for (i=0; i<bs; i++) {     /* Loop over a block of columns */
101       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
102       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
103       if (coloring->htype[0] == 'w') {
104         for (l=0; l<ncolumns[k]; l++) {
105           col            = i + bs*coloring->columns[k][l];  /* local column (in global index!) of the matrix we are probing for */
106           w3_array[col] += 1.0/dx;
107           if (i) w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */
108         }
109       } else { /* htype == 'ds' */
110         vscale_array -= cstart; /* shift pointer so global index can be used */
111         for (l=0; l<ncolumns[k]; l++) {
112           col = i + bs*coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
113           w3_array[col] += 1.0/vscale_array[col];
114           if (i) w3_array[col-1] -=  1.0/vscale_array[col-1]; /* resume original w3[col-1] */
115         }
116         vscale_array += cstart;
117       }
118       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
119       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
120 
121       /*
122        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
123                            w2 = F(x1 + dx) - F(x1)
124        */
125       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
126       ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */
127       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
128       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
129       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
130       ierr = VecResetArray(w2);CHKERRQ(ierr);
131       dy_i += nxloc; /* points to dy+i*nxloc */
132     }
133 
134     /*
135      (3-3) Loop over rows of vector, putting results into Jacobian matrix
136     */
137     nrows_k = nrows[k];
138     if (coloring->htype[0] == 'w') {
139       for (l=0; l<nrows_k; l++) {
140         row     = bs*Jentry[nz].row;   /* local row index */
141         valaddr = Jentry[nz++].valaddr;
142         spidx   = 0;
143         dy_i    = dy;
144         for (i=0; i<bs; i++) {   /* column of the block */
145           for (j=0; j<bs; j++) { /* row of the block */
146             valaddr[spidx++] = dy_i[row+j]*dx;
147           }
148           dy_i += nxloc; /* points to dy+i*nxloc */
149         }
150       }
151     } else { /* htype == 'ds' */
152       for (l=0; l<nrows_k; l++) {
153         row     = bs*Jentry[nz].row;   /* local row index */
154         col     = bs*Jentry[nz].col;   /* local column index */
155         valaddr = Jentry[nz++].valaddr;
156         spidx   = 0;
157         dy_i    = dy;
158         for (i=0; i<bs; i++) {   /* column of the block */
159           for (j=0; j<bs; j++) { /* row of the block */
160             valaddr[spidx++] = dy_i[row+j]*vscale_array[col+i];
161           }
162           dy_i += nxloc; /* points to dy+i*nxloc */
163         }
164       }
165     }
166   }
167   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169   if (vscale) {
170     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
171   }
172 
173   coloring->currentcolor = -1;
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "MatFDColoringApply_AIJ"
179 PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
180 {
181   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f;
182   PetscErrorCode ierr;
183   PetscInt       k,cstart,cend,l,row,col,nz;
184   PetscScalar    dx=0.0,*y,*xx,*w3_array;
185   PetscScalar    *vscale_array;
186   PetscReal      epsilon=coloring->error_rel,umin=coloring->umin,unorm;
187   Vec            w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale;
188   void           *fctx=coloring->fctx;
189   PetscInt       ctype=coloring->ctype,nxloc,nrows_k;
190   MatEntry       *Jentry=coloring->matentry;
191   const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows;
192 
193   PetscFunctionBegin;
194   /* create vscale for storing dx */
195   if (!vscale) {
196     if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
197       Mat_MPIAIJ     *aij=(Mat_MPIAIJ*)J->data;
198       ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr);
199     } else if (ctype == IS_COLORING_GHOSTED) {
200       ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr);
201     }
202     coloring->vscale = vscale;
203   }
204 
205   /* (1) Set w1 = F(x1) */
206   if (!coloring->fset) {
207     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
208     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
209     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
210   } else {
211     coloring->fset = PETSC_FALSE;
212   }
213 
214   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
215   if (coloring->htype[0] == 'w') {
216     /* vscale = 1./dx is a constant scalar */
217     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
218     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
219   } else {
220     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
221     ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
222     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
223     for (col=0; col<nxloc; col++) {
224       dx = xx[col];
225       if (PetscAbsScalar(dx) < umin) {
226         if (PetscRealPart(dx) >= 0.0)      dx = umin;
227         else if (PetscRealPart(dx) < 0.0 ) dx = -umin;
228       }
229       dx               *= epsilon;
230       vscale_array[col] = 1.0/dx;
231     }
232     ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
233     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
234   }
235   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
236     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
237     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
238   }
239 
240   /* (3) Loop over each color */
241   if (!coloring->w3) {
242     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
243     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
244   }
245   w3 = coloring->w3;
246 
247   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
248   if (vscale) {
249     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
250   }
251   nz   = 0;
252 
253   /*------------- reorder Jentry ----------------*/
254   PetscInt brows=coloring->brows,bcols=coloring->bcols;
255   if (bcols>1 && ctype == IS_COLORING_GHOSTED) { /* oly supported for seqaij matrix */
256     PetscInt    i,m=J->rmap->n,nbcols;
257     PetscScalar *dy=coloring->dy,*dy_k;
258 
259     printf("      use block rows impl: brows %d, bcols %d\n",brows,bcols);
260     nbcols = 0;
261     for (k=0; k<ncolors; k+=bcols) {
262       coloring->currentcolor = k;
263 
264       /*
265        (3-1) Loop over each column associated with color
266        adding the perturbation to the vector w3 = x1 + dx.
267        */
268 
269       dy_k = dy;
270       if (k + bcols > ncolors) bcols = ncolors - k;
271       for (i=0; i<bcols; i++) {
272 
273         ierr = VecCopy(x1,w3);CHKERRQ(ierr);
274         ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
275         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
276         if (coloring->htype[0] == 'w') {
277           for (l=0; l<ncolumns[k+i]; l++) {
278             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
279             w3_array[col] += 1.0/dx;
280           }
281         } else { /* htype == 'ds' */
282           vscale_array -= cstart; /* shift pointer so global index can be used */
283           for (l=0; l<ncolumns[k+i]; l++) {
284             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
285             w3_array[col] += 1.0/vscale_array[col];
286           }
287           vscale_array += cstart;
288         }
289         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
290         ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
291 
292         /*
293          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
294                            w2 = F(x1 + dx) - F(x1)
295          */
296         ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
297         ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */
298         ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
299         ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
300         ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
301         ierr = VecResetArray(w2);CHKERRQ(ierr);
302         dy_k += m; /* points to dy+i*nxloc */
303 
304       }
305 
306       /*
307        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
308        */
309       nrows_k = nrows[nbcols++];
310       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
311 
312       if (coloring->htype[0] == 'w') {
313         for (l=0; l<nrows_k; l++) {
314           row                     = Jentry[nz].row;   /* local row index */
315           *(Jentry[nz++].valaddr) = dy[row]*dx;
316         }
317       } else { /* htype == 'ds' */
318         for (l=0; l<nrows_k; l++) {
319           row                   = Jentry[nz].row;   /* local row index */
320           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
321           nz++;
322         }
323       }
324       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
325     }
326     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
328     if (vscale) {
329       ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
330     }
331 
332     coloring->currentcolor = -1;
333     PetscFunctionReturn(0);
334   } /*------------------ endof reorder Jentry ----------------*/
335 
336   for (k=0; k<ncolors; k++) {
337     coloring->currentcolor = k;
338 
339     /*
340       (3-1) Loop over each column associated with color
341       adding the perturbation to the vector w3 = x1 + dx.
342     */
343     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
344     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
345     if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
346     if (coloring->htype[0] == 'w') {
347       for (l=0; l<ncolumns[k]; l++) {
348         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
349         w3_array[col] += 1.0/dx;
350       }
351     } else { /* htype == 'ds' */
352       vscale_array -= cstart; /* shift pointer so global index can be used */
353       for (l=0; l<ncolumns[k]; l++) {
354         col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
355         w3_array[col] += 1.0/vscale_array[col];
356       }
357       vscale_array += cstart;
358     }
359     if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
360     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
361 
362     /*
363       (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
364                            w2 = F(x1 + dx) - F(x1)
365     */
366     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
367     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
368     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
369     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
370 
371     /*
372      (3-3) Loop over rows of vector, putting results into Jacobian matrix
373     */
374     nrows_k = nrows[k];
375     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
376     if (coloring->htype[0] == 'w') {
377       for (l=0; l<nrows_k; l++) {
378         row                     = Jentry[nz].row;   /* local row index */
379         *(Jentry[nz++].valaddr) = y[row]*dx;
380       }
381     } else { /* htype == 'ds' */
382       for (l=0; l<nrows_k; l++) {
383         row                   = Jentry[nz].row;   /* local row index */
384         *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
385         nz++;
386       }
387     }
388     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
389   }
390   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
391   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392   if (vscale) {
393     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
394   }
395 
396   coloring->currentcolor = -1;
397   PetscFunctionReturn(0);
398 }
399 
400 #undef __FUNCT__
401 #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ"
402 PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
403 {
404   PetscErrorCode         ierr;
405   PetscMPIInt            size,*ncolsonproc,*disp,nn;
406   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
407   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
408   PetscInt               nis=iscoloring->n,nctot,*cols;
409   IS                     *isa;
410   ISLocalToGlobalMapping map=mat->cmap->mapping;
411   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
412   Mat                    A,B;
413   PetscScalar            *A_val,*B_val,**valaddrhit;
414   MatEntry               *Jentry;
415   PetscBool              isBAIJ;
416 #if defined(PETSC_USE_CTABLE)
417   PetscTable             colmap=NULL;
418 #else
419   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
420 #endif
421 
422   PetscFunctionBegin;
423   if (ctype == IS_COLORING_GHOSTED) {
424     if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
425     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
426   }
427 
428   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
429   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
430   if (isBAIJ) {
431     Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data;
432     Mat_SeqBAIJ *spA,*spB;
433     A = aij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
434     B = aij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
435     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
436     if (!aij->colmap) {
437       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
438       colmap = aij->colmap;
439     }
440     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
441     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
442   } else {
443     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
444     Mat_SeqAIJ *spA,*spB;
445     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
446     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
447     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
448     if (!aij->colmap) {
449       /* Allow access to data structures of local part of matrix
450        - creates aij->colmap which maps global column number to local number in part B */
451       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
452       colmap = aij->colmap;
453     }
454     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
455     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
456 
457     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
458   }
459 
460   m         = mat->rmap->n/bs;
461   cstart    = mat->cmap->rstart/bs;
462   cend      = mat->cmap->rend/bs;
463   c->M      = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
464   c->N      = mat->cmap->N/bs;
465   c->m      = m;
466   c->rstart = mat->rmap->rstart/bs;
467 
468   c->ncolors = nis;
469   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
470   ierr       = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr);
471   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
472   ierr       = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
473 
474   ierr       = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr);
475   ierr       = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
476   c->matentry = Jentry;
477 
478   ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&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 = PetscMalloc(n*sizeof(PetscInt),&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,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&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 = PetscMalloc((nctot+1)*sizeof(PetscInt),&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_GHOSTED) {
517       /* Determine local number of columns of this color on this process, including ghost points */
518       nctot = n;
519       ierr  = PetscMalloc((nctot+1)*sizeof(PetscInt),&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_GHOSTED) {
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     for (j=0; j<m; j++) {
571       if (rowhit[j]) {
572         Jentry[nz].row     = j;              /* local row index */
573         Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
574         Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
575         nz++;
576       }
577     }
578     ierr = PetscFree(cols);CHKERRQ(ierr);
579   }
580   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
581 
582   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
583   if (isBAIJ) {
584     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
585     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
586     ierr = PetscMalloc(bs*mat->rmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr);
587   } else {
588     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
589     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
590   }
591 
592   if (ctype == IS_COLORING_GHOSTED) {
593     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
594   }
595   PetscFunctionReturn(0);
596 }
597