xref: /petsc/src/mat/impls/aij/mpi/fdmpiaij.c (revision ed5e4e85792117e4bf875eb7efb50a9b8564ac26)
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           /* The 'useless' ifdef is due to a bug in NVIDIA nvc 21.11, which triggers a segfault on this line. We write it in
291              another way, and it seems work. See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html
292            */
293          #if defined(PETSC_USE_COMPLEX)
294           PetscScalar *tmp = Jentry2[nz].valaddr;
295           *tmp = dy[row]*dx;
296          #else
297           *(Jentry2[nz].valaddr) = dy[row]*dx;
298          #endif
299           nz++;
300         }
301       } else { /* htype == 'ds' */
302         for (l=0; l<nrows_k; l++) {
303           row = Jentry[nz].row;   /* local row index */
304          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
305           PetscScalar *tmp = Jentry[nz].valaddr;
306           *tmp = dy[row]*vscale_array[Jentry[nz].col];
307          #else
308           *(Jentry[nz].valaddr) = dy[row]*vscale_array[Jentry[nz].col];
309          #endif
310           nz++;
311         }
312       }
313     }
314   } else { /* bcols == 1 */
315     for (k=0; k<ncolors; k++) {
316       coloring->currentcolor = k;
317 
318       /*
319        (3-1) Loop over each column associated with color
320        adding the perturbation to the vector w3 = x1 + dx.
321        */
322       ierr = VecCopy(x1,w3);CHKERRQ(ierr);
323       ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
324       if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */
325       if (coloring->htype[0] == 'w') {
326         for (l=0; l<ncolumns[k]; l++) {
327           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
328           w3_array[col] += 1.0/dx;
329         }
330       } else { /* htype == 'ds' */
331         vscale_array -= cstart; /* shift pointer so global index can be used */
332         for (l=0; l<ncolumns[k]; l++) {
333           col = coloring->columns[k][l]; /* local column (in global index!) of the matrix we are probing for */
334           w3_array[col] += 1.0/vscale_array[col];
335         }
336         vscale_array += cstart;
337       }
338       if (ctype == IS_COLORING_GLOBAL) w3_array += cstart;
339       ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
340 
341       /*
342        (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
343                            w2 = F(x1 + dx) - F(x1)
344        */
345       ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
346       ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
347       ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
348       ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
349 
350       /*
351        (3-3) Loop over rows of vector, putting results into Jacobian matrix
352        */
353       nrows_k = nrows[k];
354       ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
355       if (coloring->htype[0] == 'w') {
356         for (l=0; l<nrows_k; l++) {
357           row  = Jentry2[nz].row;   /* local row index */
358          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
359           PetscScalar *tmp = Jentry2[nz].valaddr;
360           *tmp = y[row]*dx;
361          #else
362           *(Jentry2[nz].valaddr) = y[row]*dx;
363          #endif
364           nz++;
365         }
366       } else { /* htype == 'ds' */
367         for (l=0; l<nrows_k; l++) {
368           row  = Jentry[nz].row;   /* local row index */
369          #if defined(PETSC_USE_COMPLEX) /* See https://lists.mcs.anl.gov/pipermail/petsc-users/2021-December/045158.html */
370           PetscScalar *tmp = Jentry[nz].valaddr;
371           *tmp = y[row]*vscale_array[Jentry[nz].col];
372          #else
373           *(Jentry[nz].valaddr) = y[row]*vscale_array[Jentry[nz].col];
374          #endif
375           nz++;
376         }
377       }
378       ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
379     }
380   }
381 
382 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
383   if (J->offloadmask != PETSC_OFFLOAD_UNALLOCATED) J->offloadmask = PETSC_OFFLOAD_CPU;
384 #endif
385   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
386   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
387   if (vscale) {
388     ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr);
389   }
390   coloring->currentcolor = -1;
391   if (!alreadyboundtocpu) {ierr = VecBindToCPU(x1,PETSC_FALSE);CHKERRQ(ierr);}
392   PetscFunctionReturn(0);
393 }
394 
395 PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
396 {
397   PetscErrorCode         ierr;
398   PetscMPIInt            size,*ncolsonproc,*disp,nn;
399   PetscInt               i,n,nrows,nrows_i,j,k,m,ncols,col,*rowhit,cstart,cend,colb;
400   const PetscInt         *is,*A_ci,*A_cj,*B_ci,*B_cj,*row=NULL,*ltog=NULL;
401   PetscInt               nis=iscoloring->n,nctot,*cols,tmp = 0;
402   ISLocalToGlobalMapping map=mat->cmap->mapping;
403   PetscInt               ctype=c->ctype,*spidxA,*spidxB,nz,bs,bs2,spidx;
404   Mat                    A,B;
405   PetscScalar            *A_val,*B_val,**valaddrhit;
406   MatEntry               *Jentry;
407   MatEntry2              *Jentry2;
408   PetscBool              isBAIJ,isSELL;
409   PetscInt               bcols=c->bcols;
410 #if defined(PETSC_USE_CTABLE)
411   PetscTable             colmap=NULL;
412 #else
413   PetscInt               *colmap=NULL;     /* local col number of off-diag col */
414 #endif
415 
416   PetscFunctionBegin;
417   if (ctype == IS_COLORING_LOCAL) {
418     if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping");
419     ierr = ISLocalToGlobalMappingGetIndices(map,&ltog);CHKERRQ(ierr);
420   }
421 
422   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
423   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
424   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
425   if (isBAIJ) {
426     Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)mat->data;
427     Mat_SeqBAIJ *spA,*spB;
428     A = baij->A;  spA = (Mat_SeqBAIJ*)A->data; A_val = spA->a;
429     B = baij->B;  spB = (Mat_SeqBAIJ*)B->data; B_val = spB->a;
430     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
431     if (!baij->colmap) {
432       ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
433     }
434     colmap = baij->colmap;
435     ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
436     ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
437 
438     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') {  /* create vscale for storing dx */
439       PetscInt    *garray;
440       ierr = PetscMalloc1(B->cmap->n,&garray);CHKERRQ(ierr);
441       for (i=0; i<baij->B->cmap->n/bs; i++) {
442         for (j=0; j<bs; j++) {
443           garray[i*bs+j] = bs*baij->garray[i]+j;
444         }
445       }
446       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,garray,&c->vscale);CHKERRQ(ierr);
447       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
448       ierr = PetscFree(garray);CHKERRQ(ierr);
449     }
450   } else if (isSELL) {
451     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
452     Mat_SeqSELL *spA,*spB;
453     A = sell->A;  spA = (Mat_SeqSELL*)A->data; A_val = spA->val;
454     B = sell->B;  spB = (Mat_SeqSELL*)B->data; B_val = spB->val;
455     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
456     if (!sell->colmap) {
457       /* Allow access to data structures of local part of matrix
458        - creates aij->colmap which maps global column number to local number in part B */
459       ierr = MatCreateColmap_MPISELL_Private(mat);CHKERRQ(ierr);
460     }
461     colmap = sell->colmap;
462     ierr = MatGetColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
463     ierr = MatGetColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
464 
465     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
466 
467     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
468       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,sell->garray,&c->vscale);CHKERRQ(ierr);
469       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
470     }
471   } else {
472     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
473     Mat_SeqAIJ *spA,*spB;
474     A = aij->A;  spA = (Mat_SeqAIJ*)A->data; A_val = spA->a;
475     B = aij->B;  spB = (Mat_SeqAIJ*)B->data; B_val = spB->a;
476     nz = spA->nz + spB->nz; /* total nonzero entries of mat */
477     if (!aij->colmap) {
478       /* Allow access to data structures of local part of matrix
479        - creates aij->colmap which maps global column number to local number in part B */
480       ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr);
481     }
482     colmap = aij->colmap;
483     ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
484     ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
485 
486     bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */
487 
488     if (ctype == IS_COLORING_GLOBAL && c->htype[0] == 'd') { /* create vscale for storing dx */
489       ierr = VecCreateGhost(PetscObjectComm((PetscObject)mat),mat->cmap->n,PETSC_DETERMINE,B->cmap->n,aij->garray,&c->vscale);CHKERRQ(ierr);
490       ierr = VecBindToCPU(c->vscale,PETSC_TRUE);CHKERRQ(ierr);
491     }
492   }
493 
494   m      = mat->rmap->n/bs;
495   cstart = mat->cmap->rstart/bs;
496   cend   = mat->cmap->rend/bs;
497 
498   ierr = PetscMalloc2(nis,&c->ncolumns,nis,&c->columns);CHKERRQ(ierr);
499   ierr = PetscMalloc1(nis,&c->nrows);CHKERRQ(ierr);
500   ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr);
501 
502   if (c->htype[0] == 'd') {
503     ierr        = PetscMalloc1(nz,&Jentry);CHKERRQ(ierr);
504     ierr        = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr);
505     c->matentry = Jentry;
506   } else if (c->htype[0] == 'w') {
507     ierr         = PetscMalloc1(nz,&Jentry2);CHKERRQ(ierr);
508     ierr         = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry2));CHKERRQ(ierr);
509     c->matentry2 = Jentry2;
510   } else SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"htype is not supported");
511 
512   ierr = PetscMalloc2(m+1,&rowhit,m+1,&valaddrhit);CHKERRQ(ierr);
513   nz   = 0;
514   ierr = ISColoringGetIS(iscoloring,PETSC_OWN_POINTER, PETSC_IGNORE,&c->isa);CHKERRQ(ierr);
515 
516   if (ctype == IS_COLORING_GLOBAL) {
517     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRMPI(ierr);
518     ierr = PetscMalloc2(size,&ncolsonproc,size,&disp);CHKERRQ(ierr);
519   }
520 
521   for (i=0; i<nis; i++) { /* for each local color */
522     ierr = ISGetLocalSize(c->isa[i],&n);CHKERRQ(ierr);
523     ierr = ISGetIndices(c->isa[i],&is);CHKERRQ(ierr);
524 
525     c->ncolumns[i] = n; /* local number of columns of this color on this process */
526     c->columns[i]  = (PetscInt*)is;
527 
528     if (ctype == IS_COLORING_GLOBAL) {
529       /* Determine nctot, the total (parallel) number of columns of this color */
530       /* ncolsonproc[j]: local ncolumns on proc[j] of this color */
531       ierr  = PetscMPIIntCast(n,&nn);CHKERRQ(ierr);
532       ierr  = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
533       nctot = 0; for (j=0; j<size; j++) nctot += ncolsonproc[j];
534       if (!nctot) {
535         ierr = PetscInfo(mat,"Coloring of matrix has some unneeded colors with no corresponding rows\n");CHKERRQ(ierr);
536       }
537 
538       disp[0] = 0;
539       for (j=1; j<size; j++) {
540         disp[j] = disp[j-1] + ncolsonproc[j-1];
541       }
542 
543       /* Get cols, the complete list of columns for this color on each process */
544       ierr = PetscMalloc1(nctot+1,&cols);CHKERRQ(ierr);
545       ierr = MPI_Allgatherv((void*)is,n,MPIU_INT,cols,ncolsonproc,disp,MPIU_INT,PetscObjectComm((PetscObject)mat));CHKERRMPI(ierr);
546     } else if (ctype == IS_COLORING_LOCAL) {
547       /* Determine local number of columns of this color on this process, including ghost points */
548       nctot = n;
549       cols  = (PetscInt*)is;
550     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not provided for this MatFDColoring type");
551 
552     /* Mark all rows affect by these columns */
553     ierr    = PetscArrayzero(rowhit,m);CHKERRQ(ierr);
554     bs2     = bs*bs;
555     nrows_i = 0;
556     for (j=0; j<nctot; j++) { /* loop over columns*/
557       if (ctype == IS_COLORING_LOCAL) {
558         col = ltog[cols[j]];
559       } else {
560         col = cols[j];
561       }
562       if (col >= cstart && col < cend) { /* column is in A, diagonal block of mat */
563         tmp      = A_ci[col-cstart];
564         row      = A_cj + tmp;
565         nrows    = A_ci[col-cstart+1] - tmp;
566         nrows_i += nrows;
567 
568         /* loop over columns of A marking them in rowhit */
569         for (k=0; k<nrows; k++) {
570           /* set valaddrhit for part A */
571           spidx            = bs2*spidxA[tmp + k];
572           valaddrhit[*row] = &A_val[spidx];
573           rowhit[*row++]   = col - cstart + 1; /* local column index */
574         }
575       } else { /* column is in B, off-diagonal block of mat */
576 #if defined(PETSC_USE_CTABLE)
577         ierr = PetscTableFind(colmap,col+1,&colb);CHKERRQ(ierr);
578         colb--;
579 #else
580         colb = colmap[col] - 1; /* local column index */
581 #endif
582         if (colb == -1) {
583           nrows = 0;
584         } else {
585           colb  = colb/bs;
586           tmp   = B_ci[colb];
587           row   = B_cj + tmp;
588           nrows = B_ci[colb+1] - tmp;
589         }
590         nrows_i += nrows;
591         /* loop over columns of B marking them in rowhit */
592         for (k=0; k<nrows; k++) {
593           /* set valaddrhit for part B */
594           spidx            = bs2*spidxB[tmp + k];
595           valaddrhit[*row] = &B_val[spidx];
596           rowhit[*row++]   = colb + 1 + cend - cstart; /* local column index */
597         }
598       }
599     }
600     c->nrows[i] = nrows_i;
601 
602     if (c->htype[0] == 'd') {
603       for (j=0; j<m; j++) {
604         if (rowhit[j]) {
605           Jentry[nz].row     = j;              /* local row index */
606           Jentry[nz].col     = rowhit[j] - 1;  /* local column index */
607           Jentry[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
608           nz++;
609         }
610       }
611     } else { /* c->htype == 'wp' */
612       for (j=0; j<m; j++) {
613         if (rowhit[j]) {
614           Jentry2[nz].row     = j;              /* local row index */
615           Jentry2[nz].valaddr = valaddrhit[j];  /* address of mat value for this entry */
616           nz++;
617         }
618       }
619     }
620     if (ctype == IS_COLORING_GLOBAL) {
621       ierr = PetscFree(cols);CHKERRQ(ierr);
622     }
623   }
624   if (ctype == IS_COLORING_GLOBAL) {
625     ierr = PetscFree2(ncolsonproc,disp);CHKERRQ(ierr);
626   }
627 
628   if (bcols > 1) { /* reorder Jentry for faster MatFDColoringApply() */
629     ierr = MatFDColoringSetUpBlocked_AIJ_Private(mat,c,nz);CHKERRQ(ierr);
630   }
631 
632   if (isBAIJ) {
633     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
634     ierr = MatRestoreColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
635     ierr = PetscMalloc1(bs*mat->rmap->n,&c->dy);CHKERRQ(ierr);
636   } else if (isSELL) {
637     ierr = MatRestoreColumnIJ_SeqSELL_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
638     ierr = MatRestoreColumnIJ_SeqSELL_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
639   } else {
640     ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr);
641     ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr);
642   }
643 
644   ierr = ISColoringRestoreIS(iscoloring,PETSC_OWN_POINTER,&c->isa);CHKERRQ(ierr);
645   ierr = PetscFree2(rowhit,valaddrhit);CHKERRQ(ierr);
646 
647   if (ctype == IS_COLORING_LOCAL) {
648     ierr = ISLocalToGlobalMappingRestoreIndices(map,&ltog);CHKERRQ(ierr);
649   }
650   ierr = PetscInfo3(c,"ncolors %D, brows %D and bcols %D are used.\n",c->ncolors,c->brows,c->bcols);CHKERRQ(ierr);
651   PetscFunctionReturn(0);
652 }
653 
654 PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c)
655 {
656   PetscErrorCode ierr;
657   PetscInt       bs,nis=iscoloring->n,m=mat->rmap->n;
658   PetscBool      isBAIJ,isSELL;
659 
660   PetscFunctionBegin;
661   /* set default brows and bcols for speedup inserting the dense matrix into sparse Jacobian;
662    bcols is chosen s.t. dy-array takes 50% of memory space as mat */
663   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
664   ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);
665   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPISELL,&isSELL);CHKERRQ(ierr);
666   if (isBAIJ || m == 0) {
667     c->brows = m;
668     c->bcols = 1;
669   } else if (isSELL) {
670     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
671     Mat_MPISELL *sell=(Mat_MPISELL*)mat->data;
672     Mat_SeqSELL *spA,*spB;
673     Mat        A,B;
674     PetscInt   nz,brows,bcols;
675     PetscReal  mem;
676 
677     bs    = 1; /* only bs=1 is supported for MPISELL matrix */
678 
679     A = sell->A;  spA = (Mat_SeqSELL*)A->data;
680     B = sell->B;  spB = (Mat_SeqSELL*)B->data;
681     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
682     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
683     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
684     brows = 1000/bcols;
685     if (bcols > nis) bcols = nis;
686     if (brows == 0 || brows > m) brows = m;
687     c->brows = brows;
688     c->bcols = bcols;
689   } else { /* mpiaij matrix */
690     /* bcols is chosen s.t. dy-array takes 50% of local memory space as mat */
691     Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data;
692     Mat_SeqAIJ *spA,*spB;
693     Mat        A,B;
694     PetscInt   nz,brows,bcols;
695     PetscReal  mem;
696 
697     bs    = 1; /* only bs=1 is supported for MPIAIJ matrix */
698 
699     A = aij->A;  spA = (Mat_SeqAIJ*)A->data;
700     B = aij->B;  spB = (Mat_SeqAIJ*)B->data;
701     nz = spA->nz + spB->nz; /* total local nonzero entries of mat */
702     mem = nz*(sizeof(PetscScalar) + sizeof(PetscInt)) + 3*m*sizeof(PetscInt);
703     bcols = (PetscInt)(0.5*mem /(m*sizeof(PetscScalar)));
704     brows = 1000/bcols;
705     if (bcols > nis) bcols = nis;
706     if (brows == 0 || brows > m) brows = m;
707     c->brows = brows;
708     c->bcols = bcols;
709   }
710 
711   c->M       = mat->rmap->N/bs;         /* set the global rows and columns and local rows */
712   c->N       = mat->cmap->N/bs;
713   c->m       = mat->rmap->n/bs;
714   c->rstart  = mat->rmap->rstart/bs;
715   c->ncolors = nis;
716   PetscFunctionReturn(0);
717 }
718 
719 /*@C
720 
721     MatFDColoringSetValues - takes a matrix in compressed color format and enters the matrix into a PETSc Mat
722 
723    Collective on J
724 
725    Input Parameters:
726 +    J - the sparse matrix
727 .    coloring - created with MatFDColoringCreate() and a local coloring
728 -    y - column major storage of matrix values with one color of values per column, the number of rows of y should match
729          the number of local rows of J and the number of columns is the number of colors.
730 
731    Level: intermediate
732 
733    Notes: the matrix in compressed color format may come from an Automatic Differentiation code
734 
735    The code will be slightly faster if MatFDColoringSetBlockSize(coloring,PETSC_DEFAULT,nc); is called immediately after creating the coloring
736 
737 .seealso: MatFDColoringCreate(), ISColoring, ISColoringCreate(), ISColoringSetType(), IS_COLORING_LOCAL, MatFDColoringSetBlockSize()
738 
739 @*/
740 PetscErrorCode  MatFDColoringSetValues(Mat J,MatFDColoring coloring,const PetscScalar *y)
741 {
742   PetscErrorCode    ierr;
743   MatEntry2         *Jentry2;
744   PetscInt          row,i,nrows_k,l,ncolors,nz = 0,bcols,nbcols = 0;
745   const PetscInt    *nrows;
746   PetscBool         eq;
747 
748   PetscFunctionBegin;
749   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
750   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
751   ierr = PetscObjectCompareId((PetscObject)J,coloring->matid,&eq);CHKERRQ(ierr);
752   if (!eq) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Matrix used with MatFDColoringSetValues() must be that used with MatFDColoringCreate()");
753   Jentry2 = coloring->matentry2;
754   nrows   = coloring->nrows;
755   ncolors = coloring->ncolors;
756   bcols   = coloring->bcols;
757 
758   for (i=0; i<ncolors; i+=bcols) {
759     nrows_k = nrows[nbcols++];
760     for (l=0; l<nrows_k; l++) {
761       row                      = Jentry2[nz].row;   /* local row index */
762       *(Jentry2[nz++].valaddr) = y[row];
763     }
764     y += bcols*coloring->m;
765   }
766   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
768   PetscFunctionReturn(0);
769 }
770