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