xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
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 = VecBindToCPU(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 intensively in particular piece of scalar CPU code; won't benefit from bouncing back and forth to the GPU */
65     ierr = VecBindToCPU(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 = VecBindToCPU(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   PetscBool         alreadyboundtocpu;
180 
181   PetscFunctionBegin;
182   ierr = VecBoundToCPU(x1,&alreadyboundtocpu);CHKERRQ(ierr);
183   ierr = VecBindToCPU(x1,PETSC_TRUE);CHKERRQ(ierr);
184   if ((ctype == IS_COLORING_LOCAL) && (J->ops->fdcoloringapply == MatFDColoringApply_AIJ)) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_SUP,"Must call MatColoringUseDM() with IS_COLORING_LOCAL");
185   /* (1) Set w1 = F(x1) */
186   if (!coloring->fset) {
187     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
188     ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
189     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
190   } else {
191     coloring->fset = PETSC_FALSE;
192   }
193 
194   /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */
195   if (coloring->htype[0] == 'w') {
196     /* vscale = 1./dx is a constant scalar */
197     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
198     dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon);
199   } else {
200     ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr);
201     ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
202     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
203     for (col=0; col<nxloc; col++) {
204       dx = xx[col];
205       if (PetscAbsScalar(dx) < umin) {
206         if (PetscRealPart(dx) >= 0.0)      dx = umin;
207         else if (PetscRealPart(dx) < 0.0) dx = -umin;
208       }
209       dx               *= epsilon;
210       vscale_array[col] = 1.0/dx;
211     }
212     ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
213     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
214   }
215   if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') {
216     ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
217     ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
218   }
219 
220   /* (3) Loop over each color */
221   if (!coloring->w3) {
222     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
223     ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr);
224   }
225   w3 = coloring->w3;
226 
227   ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */
228   if (vscale) {
229     ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr);
230   }
231   nz = 0;
232 
233   if (coloring->bcols > 1) { /* use blocked insertion of Jentry */
234     PetscInt    i,m=J->rmap->n,nbcols,bcols=coloring->bcols;
235     PetscScalar *dy=coloring->dy,*dy_k;
236 
237     nbcols = 0;
238     for (k=0; k<ncolors; k+=bcols) {
239 
240       /*
241        (3-1) Loop over each column associated with color
242        adding the perturbation to the vector w3 = x1 + dx.
243        */
244 
245       dy_k = dy;
246       if (k + bcols > ncolors) bcols = ncolors - k;
247       for (i=0; i<bcols; i++) {
248         coloring->currentcolor = k+i;
249 
250         ierr = VecCopy(x1,w3);CHKERRQ(ierr);
251         ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
252         if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
253         if (coloring->htype[0] == 'w') {
254           for (l=0; l<ncolumns[k+i]; l++) {
255             col = coloring->columns[k+i][l]; /* local column (in global index!) of the matrix we are probing for */
256             w3_array[col] += 1.0/dx;
257           }
258         } else { /* htype == 'ds' */
259           vscale_array -= cstart; /* shift pointer so global index can be used */
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/vscale_array[col];
263           }
264           vscale_array += cstart;
265         }
266         if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
267         ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
268 
269         /*
270          (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
271                            w2 = F(x1 + dx) - F(x1)
272          */
273         ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
274         ierr = VecPlaceArray(w2,dy_k);CHKERRQ(ierr); /* place w2 to the array dy_i */
275         ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
276         ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
277         ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
278         ierr = VecResetArray(w2);CHKERRQ(ierr);
279         dy_k += m; /* points to dy+i*nxloc */
280       }
281 
282       /*
283        (3-3) Loop over block rows of vector, putting results into Jacobian matrix
284        */
285       nrows_k = nrows[nbcols++];
286 
287       if (coloring->htype[0] == 'w') {
288         for (l=0; l<nrows_k; l++) {
289           row                      = Jentry2[nz].row;   /* local row index */
290           *(Jentry2[nz++].valaddr) = dy[row]*dx;
291         }
292       } else { /* htype == 'ds' */
293         for (l=0; l<nrows_k; l++) {
294           row                   = Jentry[nz].row;   /* local row index */
295           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
296           nz++;
297         }
298       }
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 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
358   if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU;
359 #endif
360   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
361   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
362   if (vscale) {
363     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
364   }
365   coloring->currentcolor = -1;
366   if (!alreadyboundtocpu) {ierr = VecBindToCPU(x1,PETSC_FALSE);CHKERRQ(ierr);}
367   PetscFunctionReturn(0);
368 }
369 
370 PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
371 {
372   PetscErrorCode         ierr;
373   PetscMPIInt            size,*ncolsonproc,*disp,nn;
374   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
375   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
376   PetscInt               nis=iscoloring->n,nctot,*cols,tmp = 0;
377   ISLocalToGlobalMapping map=mat->cmap->mapping;
378   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
379   Mat                    A,B;
380   PetscScalar            *A_val,*B_val,**valaddrhit;
381   MatEntry               *Jentry;
382   MatEntry2              *Jentry2;
383   PetscBool              isBAIJ,isSELL;
384   PetscInt               bcols=c->bcols;
385 #if defined(PETSC_USE_CTABLE)
386   PetscTable             colmap=NULL;
387 #else
388   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
389 #endif
390 
391   PetscFunctionBegin;
392   if (ctype == IS_COLORING_LOCAL) {
393     if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
394     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
395   }
396 
397   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
398   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
399   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
400   if (isBAIJ) {
401     Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data;
402     Mat_SeqBAIJ *spA,*spB;
403     A = baij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
404     B = baij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
405     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
406     if (!baij->colmap) {
407       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
408     }
409     colmap = baij->colmap;
410     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
411     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
412 
413     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') {  /* create vscale for storing dx */
414       PetscInt    *garray;
415       ierr = PetscMalloc1(B->cmap->n,&garray);CHKERRQ(ierr);
416       for (i=0; i<baij->B->cmap->n/bs; i++) {
417         for (j=0; j<bs; j++) {
418           garray[i*bs+j] = bs*baij->garray[i]+j;
419         }
420       }
421       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr);
422       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
423       ierr = PetscFree(garray);CHKERRQ(ierr);
424     }
425   } else if (isSELL) {
426     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
427     Mat_SeqSELL *spA,*spB;
428     A = sell->A;  spA = (Mat_SeqSELL*)A->data; A_val = spA->val;
429     B = sell->B;  spB = (Mat_SeqSELL*)B->data; B_val = spB->val;
430     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
431     if (!sell->colmap) {
432       /* Allow access to data structures of local part of matrix
433        - creates aij->colmap which maps global column number to local number in part B */
434       ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr);
435     }
436     colmap = sell->colmap;
437     ierr = MatGetColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
438     ierr = MatGetColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
439 
440     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
441 
442     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
443       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,sell->garray,&c->vscale);CHKERRQ(ierr);
444       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
445     }
446   } else {
447     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
448     Mat_SeqAIJ *spA,*spB;
449     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
450     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
451     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
452     if (!aij->colmap) {
453       /* Allow access to data structures of local part of matrix
454        - creates aij->colmap which maps global column number to local number in part B */
455       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
456     }
457     colmap = aij->colmap;
458     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
459     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
460 
461     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
462 
463     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
464       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
465       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
466     }
467   }
468 
469   m      = mat->rmap->n/bs;
470   cstart = mat->cmap->rstart/bs;
471   cend   = mat->cmap->rend/bs;
472 
473   ierr = PetscMalloc2(nis,&c->ncolumns,nis,&c->columns);CHKERRQ(ierr);
474   ierr = PetscMalloc1(nis,&c->nrows);CHKERRQ(ierr);
475   ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
476 
477   if (c->htype[0] == 'd') {
478     ierr        = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr);
479     ierr        = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
480     c->matentry = Jentry;
481   } else if (c->htype[0] == 'w') {
482     ierr         = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr);
483     ierr         = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr);
484     c->matentry2 = Jentry2;
485   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
486 
487   ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr);
488   nz   = 0;
489   ierr = ISColoringGetIS(iscoloring,PETSC_OWN_POINTER, PETSC_IGNORE,&c->isa);CHKERRQ(ierr);
490 
491   if (ctype == IS_COLORING_GLOBAL) {
492     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRMPI(ierr);
493     ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr);
494   }
495 
496   for (i=0; i<nis; i++) { /* for each local color */
497     ierr = ISGetLocalSize(c->isa[i],&n);CHKERRQ(ierr);
498     ierr = ISGetIndices(c->isa[i],&is);CHKERRQ(ierr);
499 
500     c->ncolumns[i] = n; /* local number of columns of this color on this process */
501     c->columns[i]  = (PetscInt*)is;
502 
503     if (ctype == IS_COLORING_GLOBAL) {
504       /* Determine nctot, the total (parallel) number of columns of this color */
505       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
506       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
507       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
508       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
509       if (!nctot) {
510         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
511       }
512 
513       disp[0] = 0;
514       for (j=1; j<size; j++) {
515         disp[j] = disp[j-1] + ncolsonproc[j-1];
516       }
517 
518       /* Get cols, the complete list of columns for this color on each process */
519       ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr);
520       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
521     } else if (ctype == IS_COLORING_LOCAL) {
522       /* Determine local number of columns of this color on this process, including ghost points */
523       nctot = n;
524       cols  = (PetscInt*)is;
525     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
526 
527     /* Mark all rows affect by these columns */
528     ierr    = PetscArrayzero(rowhit,m);CHKERRQ(ierr);
529     bs2     = bs*bs;
530     nrows_i = 0;
531     for (j=0; j<nctot; j++) { /* loop over columns*/
532       if (ctype == IS_COLORING_LOCAL) {
533         col = ltog[cols[j]];
534       } else {
535         col = cols[j];
536       }
537       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
538         tmp      = A_ci[col-cstart];
539         row      = A_cj + tmp;
540         nrows    = A_ci[col-cstart+1] - tmp;
541         nrows_i += nrows;
542 
543         /* loop over columns of A marking them in rowhit */
544         for (k=0; k<nrows; k++) {
545           /* set valaddrhit for part A */
546           spidx            = bs2*spidxA[tmp + k];
547           valaddrhit[*row] = &A_val[spidx];
548           rowhit[*row++]   = col - cstart + 1; /* local column index */
549         }
550       } else { /* column is in B, off-diagonal block of mat */
551 #if defined(PETSC_USE_CTABLE)
552         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
553         colb--;
554 #else
555         colb = colmap[col] - 1; /* local column index */
556 #endif
557         if (colb == -1) {
558           nrows = 0;
559         } else {
560           colb  = colb/bs;
561           tmp   = B_ci[colb];
562           row   = B_cj + tmp;
563           nrows = B_ci[colb+1] - tmp;
564         }
565         nrows_i += nrows;
566         /* loop over columns of B marking them in rowhit */
567         for (k=0; k<nrows; k++) {
568           /* set valaddrhit for part B */
569           spidx            = bs2*spidxB[tmp + k];
570           valaddrhit[*row] = &B_val[spidx];
571           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
572         }
573       }
574     }
575     c->nrows[i] = nrows_i;
576 
577     if (c->htype[0] == 'd') {
578       for (j=0; j<m; j++) {
579         if (rowhit[j]) {
580           Jentry[nz].row     = j;              /* local row index */
581           Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
582           Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
583           nz++;
584         }
585       }
586     } else { /* c->htype == 'wp' */
587       for (j=0; j<m; j++) {
588         if (rowhit[j]) {
589           Jentry2[nz].row     = j;              /* local row index */
590           Jentry2[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
591           nz++;
592         }
593       }
594     }
595     if (ctype == IS_COLORING_GLOBAL) {
596       ierr = PetscFree(cols);CHKERRQ(ierr);
597     }
598   }
599   if (ctype == IS_COLORING_GLOBAL) {
600     ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
601   }
602 
603   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
604     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
605   }
606 
607   if (isBAIJ) {
608     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
609     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
610     ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr);
611   } else if (isSELL) {
612     ierr = MatRestoreColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
613     ierr = MatRestoreColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
614   } else {
615     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
616     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
617   }
618 
619   ierr = ISColoringRestoreIS(iscoloring,PETSC_OWN_POINTER,&c->isa);CHKERRQ(ierr);
620   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
621 
622   if (ctype == IS_COLORING_LOCAL) {
623     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
624   }
625   ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
626   PetscFunctionReturn(0);
627 }
628 
629 PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
630 {
631   PetscErrorCode ierr;
632   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
633   PetscBool      isBAIJ,isSELL;
634 
635   PetscFunctionBegin;
636   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
637    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
638   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
639   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
640   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
641   if (isBAIJ || m == 0) {
642     c->brows = m;
643     c->bcols = 1;
644   } else if (isSELL) {
645     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
646     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
647     Mat_SeqSELL *spA,*spB;
648     Mat        A,B;
649     PetscInt   nz,brows,bcols;
650     PetscReal  mem;
651 
652     bs    = 1; /* only bs=1 is supported for MPISELL matrix */
653 
654     A = sell->A;  spA = (Mat_SeqSELL*)A->data;
655     B = sell->B;  spB = (Mat_SeqSELL*)B->data;
656     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
657     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
658     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
659     brows = 1000/bcols;
660     if (bcols > nis) bcols = nis;
661     if (brows == 0 || brows > m) brows = m;
662     c->brows = brows;
663     c->bcols = bcols;
664   } else { /* mpiaij matrix */
665     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
666     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
667     Mat_SeqAIJ *spA,*spB;
668     Mat        A,B;
669     PetscInt   nz,brows,bcols;
670     PetscReal  mem;
671 
672     bs    = 1; /* only bs=1 is supported for MPIAIJ matrix */
673 
674     A = aij->A;  spA = (Mat_SeqAIJ*)A->data;
675     B = aij->B;  spB = (Mat_SeqAIJ*)B->data;
676     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
677     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
678     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
679     brows = 1000/bcols;
680     if (bcols > nis) bcols = nis;
681     if (brows == 0 || brows > m) brows = m;
682     c->brows = brows;
683     c->bcols = bcols;
684   }
685 
686   c->M       = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
687   c->N       = mat->cmap->N/bs;
688   c->m       = mat->rmap->n/bs;
689   c->rstart  = mat->rmap->rstart/bs;
690   c->ncolors = nis;
691   PetscFunctionReturn(0);
692 }
693 
694 /*@C
695 
696     MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc Mat
697 
698    Collective on J
699 
700    Input Parameters:
701 +    J - the sparse matrix
702 .    coloring - created with MatFDColoringCreate() and a local coloring
703 -    y - column major storage of matrix values with one color of values per column, the number of rows of y should match
704          the number of local rows of J and the number of columns is the number of colors.
705 
706    Level: intermediate
707 
708    Notes: the matrix in compressed color format may come from an Automatic Differentiation code
709 
710    The code will be slightly faster if MatFDColoringSetBlockSize(coloring,PETSC_DEFAULT,nc); is called immediately after creating the coloring
711 
712 .seealso: MatFDColoringCreate(), ISColoring, ISColoringCreate(), ISColoringSetType(), IS_COLORING_LOCAL, MatFDColoringSetBlockSize()
713 
714 @*/
715 PetscErrorCode  MatFDColoringSetValues(Mat J,MatFDColoring coloring,const PetscScalar *y)
716 {
717   PetscErrorCode    ierr;
718   MatEntry2         *Jentry2;
719   PetscInt          row,i,nrows_k,l,ncolors,nz = 0,bcols,nbcols = 0;
720   const PetscInt    *nrows;
721   PetscBool         eq;
722 
723   PetscFunctionBegin;
724   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
725   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
726   ierr = PetscObjectCompareId((PetscObject)J,coloring->matid,&eq);CHKERRQ(ierr);
727   if (!eq) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()");
728   Jentry2 = coloring->matentry2;
729   nrows   = coloring->nrows;
730   ncolors = coloring->ncolors;
731   bcols   = coloring->bcols;
732 
733   for (i=0; i<ncolors; i+=bcols) {
734     nrows_k = nrows[nbcols++];
735     for (l=0; l<nrows_k; l++) {
736       row                      = Jentry2[nz].row;   /* local row index */
737       *(Jentry2[nz++].valaddr) = y[row];
738     }
739     y += bcols*coloring->m;
740   }
741   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
742   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
743   PetscFunctionReturn(0);
744 }
745