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