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