xref: /petsc/src/ts/tutorials/autodiff/adolc-utils/sparse.cxx (revision 2da392cc7c10228af19ad9843ce5155178acb644)
1 #include <petscmat.h>
2 #include <adolc/adolc_sparse.h>
3 
4 /*
5    REQUIRES configuration of PETSc with option --download-adolc.
6 
7    For documentation on ADOL-C, see
8      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
9 */
10 
11 /*
12   Basic printing for sparsity pattern
13 
14   Input parameters:
15   comm     - MPI communicator
16   m        - number of rows
17 
18   Output parameter:
19   sparsity - matrix sparsity pattern, typically computed using an ADOL-C function such as jac_pat
20 */
21 PetscErrorCode PrintSparsity(MPI_Comm comm,PetscInt m,unsigned int **sparsity)
22 {
23   PetscErrorCode ierr;
24   PetscInt       i,j;
25 
26   PetscFunctionBegin;
27   ierr = PetscPrintf(comm,"Sparsity pattern:\n");CHKERRQ(ierr);
28   for (i=0; i<m ;i++) {
29     ierr = PetscPrintf(comm,"\n %2d: ",i);CHKERRQ(ierr);
30     for (j=1; j<= (PetscInt) sparsity[i][0] ;j++)
31       ierr = PetscPrintf(comm," %2d ",sparsity[i][j]);CHKERRQ(ierr);
32   }
33   ierr = PetscPrintf(comm,"\n\n");CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 /*
38   Generate a seed matrix defining the partition of columns of a matrix by a particular coloring,
39   used for matrix compression
40 
41   Input parameter:
42   iscoloring - the index set coloring to be used
43 
44   Output parameter:
45   S          - the resulting seed matrix
46 
47   Notes:
48   Before calling GenerateSeedMatrix, Seed should be allocated as a logically 2d array with number of
49   rows equal to the matrix to be compressed and number of columns equal to the number of colors used
50   in iscoloring.
51 */
52 PetscErrorCode GenerateSeedMatrix(ISColoring iscoloring,PetscScalar **S)
53 {
54   PetscErrorCode ierr;
55   IS             *is;
56   PetscInt       p,size,colour,j;
57   const PetscInt *indices;
58 
59   PetscFunctionBegin;
60   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&p,&is);CHKERRQ(ierr);
61   for (colour=0; colour<p; colour++) {
62     ierr = ISGetLocalSize(is[colour],&size);CHKERRQ(ierr);
63     ierr = ISGetIndices(is[colour],&indices);CHKERRQ(ierr);
64     for (j=0; j<size; j++)
65       S[indices[j]][colour] = 1.;
66     ierr = ISRestoreIndices(is[colour],&indices);CHKERRQ(ierr);
67   }
68   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 /*
73   Establish a look-up vector whose entries contain the colour used for that diagonal entry. Clearly
74   we require the number of dependent and independent variables to be equal in this case.
75   Input parameters:
76   S        - the seed matrix defining the coloring
77   sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
78              function, such as jac_pat or hess_pat
79   m        - the number of rows of Seed (and the matrix to be recovered)
80   p        - the number of colors used (also the number of columns in Seed)
81   Output parameter:
82   R        - the recovery vector to be used for de-compression
83 */
84 PetscErrorCode GenerateSeedMatrixPlusRecovery(ISColoring iscoloring,PetscScalar **S,PetscScalar *R)
85 {
86   PetscErrorCode ierr;
87   IS             *is;
88   PetscInt       p,size,colour,j;
89   const PetscInt *indices;
90 
91   PetscFunctionBegin;
92   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&p,&is);CHKERRQ(ierr);
93   for (colour=0; colour<p; colour++) {
94     ierr = ISGetLocalSize(is[colour],&size);CHKERRQ(ierr);
95     ierr = ISGetIndices(is[colour],&indices);CHKERRQ(ierr);
96     for (j=0; j<size; j++) {
97       S[indices[j]][colour] = 1.;
98       R[indices[j]] = colour;
99     }
100     ierr = ISRestoreIndices(is[colour],&indices);CHKERRQ(ierr);
101   }
102   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 /*
107   Establish a look-up matrix whose entries contain the column coordinates of the corresponding entry
108   in a matrix which has been compressed using the coloring defined by some seed matrix
109 
110   Input parameters:
111   S        - the seed matrix defining the coloring
112   sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
113              function, such as jac_pat or hess_pat
114   m        - the number of rows of Seed (and the matrix to be recovered)
115   p        - the number of colors used (also the number of columns in Seed)
116 
117   Output parameter:
118   R        - the recovery matrix to be used for de-compression
119 */
120 PetscErrorCode GetRecoveryMatrix(PetscScalar **S,unsigned int **sparsity,PetscInt m,PetscInt p,PetscScalar **R)
121 {
122   PetscInt i,j,k,colour;
123 
124   PetscFunctionBegin;
125   for (i=0; i<m; i++) {
126     for (colour=0; colour<p; colour++) {
127       R[i][colour] = -1.;
128       for (k=1; k<=(PetscInt) sparsity[i][0]; k++) {
129         j = (PetscInt) sparsity[i][k];
130         if (S[j][colour] == 1.) {
131           R[i][colour] = j;
132           break;
133         }
134       }
135     }
136   }
137   PetscFunctionReturn(0);
138 }
139 
140 /*
141   Recover the values of a sparse matrix from a compressed format and insert these into a matrix
142 
143   Input parameters:
144   mode - use INSERT_VALUES or ADD_VALUES, as required
145   m    - number of rows of matrix.
146   p    - number of colors used in the compression of J (also the number of columns of R)
147   R    - recovery matrix to use in the decompression procedure
148   C    - compressed matrix to recover values from
149   a    - shift value for implicit problems (select NULL or unity for explicit problems)
150 
151   Output parameter:
152   A    - Mat to be populated with values from compressed matrix
153 */
154 PetscErrorCode RecoverJacobian(Mat A,InsertMode mode,PetscInt m,PetscInt p,PetscScalar **R,PetscScalar **C,PetscReal *a)
155 {
156   PetscErrorCode ierr;
157   PetscInt       i,j,colour;
158 
159   PetscFunctionBegin;
160   for (i=0; i<m; i++) {
161     for (colour=0; colour<p; colour++) {
162       j = (PetscInt) R[i][colour];
163       if (j != -1) {
164         if (a)
165           C[i][colour] *= *a;
166         ierr = MatSetValues(A,1,&i,1,&j,&C[i][colour],mode);CHKERRQ(ierr);
167       }
168     }
169   }
170   PetscFunctionReturn(0);
171 }
172 
173 /*
174   Recover the values of the local portion of a sparse matrix from a compressed format and insert
175   these into the local portion of a matrix
176 
177   Input parameters:
178   mode - use INSERT_VALUES or ADD_VALUES, as required
179   m    - number of rows of matrix.
180   p    - number of colors used in the compression of J (also the number of columns of R)
181   R    - recovery matrix to use in the decompression procedure
182   C    - compressed matrix to recover values from
183   a    - shift value for implicit problems (select NULL or unity for explicit problems)
184 
185   Output parameter:
186   A    - Mat to be populated with values from compressed matrix
187 */
188 PetscErrorCode RecoverJacobianLocal(Mat A,InsertMode mode,PetscInt m,PetscInt p,PetscScalar **R,PetscScalar **C,PetscReal *a)
189 {
190   PetscErrorCode ierr;
191   PetscInt       i,j,colour;
192 
193   PetscFunctionBegin;
194   for (i=0; i<m; i++) {
195     for (colour=0; colour<p; colour++) {
196       j = (PetscInt) R[i][colour];
197       if (j != -1) {
198         if (a)
199           C[i][colour] *= *a;
200         ierr = MatSetValuesLocal(A,1,&i,1,&j,&C[i][colour],mode);CHKERRQ(ierr);
201       }
202     }
203   }
204   PetscFunctionReturn(0);
205 }
206 
207 /*
208   Recover the diagonal of the Jacobian from its compressed matrix format
209 
210   Input parameters:
211   mode - use INSERT_VALUES or ADD_VALUES, as required
212   m    - number of rows of matrix.
213   R    - recovery vector to use in the decompression procedure
214   C    - compressed matrix to recover values from
215   a    - shift value for implicit problems (select NULL or unity for explicit problems)
216 
217   Output parameter:
218   diag - Vec to be populated with values from compressed matrix
219 */
220 PetscErrorCode RecoverDiagonal(Vec diag,InsertMode mode,PetscInt m,PetscScalar *R,PetscScalar **C,PetscReal *a)
221 {
222   PetscErrorCode ierr;
223   PetscInt       i,colour;
224 
225   PetscFunctionBegin;
226   for (i=0; i<m; i++) {
227     colour = (PetscInt)R[i];
228     if (a)
229       C[i][colour] *= *a;
230     ierr = VecSetValues(diag,1,&i,&C[i][colour],mode);CHKERRQ(ierr);
231   }
232   PetscFunctionReturn(0);
233 }
234 
235 /*
236   Recover the local portion of the diagonal of the Jacobian from its compressed matrix format
237 
238   Input parameters:
239   mode - use INSERT_VALUES or ADD_VALUES, as required
240   m    - number of rows of matrix.
241   R    - recovery vector to use in the decompression procedure
242   C    - compressed matrix to recover values from
243   a    - shift value for implicit problems (select NULL or unity for explicit problems)
244 
245   Output parameter:
246   diag - Vec to be populated with values from compressed matrix
247 */
248 PetscErrorCode RecoverDiagonalLocal(Vec diag,InsertMode mode,PetscInt m,PetscScalar *R,PetscScalar **C,PetscReal *a)
249 {
250   PetscErrorCode ierr;
251   PetscInt       i,colour;
252 
253   PetscFunctionBegin;
254   for (i=0; i<m; i++) {
255     colour = (PetscInt)R[i];
256     if (a)
257       C[i][colour] *= *a;
258     ierr = VecSetValuesLocal(diag,1,&i,&C[i][colour],mode);CHKERRQ(ierr);
259   }
260   PetscFunctionReturn(0);
261 }
262