1c4762a1bSJed Brown #include <petscmat.h>
2c4762a1bSJed Brown #include <adolc/adolc_sparse.h>
3c4762a1bSJed Brown
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc.
6c4762a1bSJed Brown
7c4762a1bSJed Brown For documentation on ADOL-C, see
8c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
9c4762a1bSJed Brown */
10c4762a1bSJed Brown
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown Basic printing for sparsity pattern
13c4762a1bSJed Brown
14c4762a1bSJed Brown Input parameters:
15c4762a1bSJed Brown comm - MPI communicator
16c4762a1bSJed Brown m - number of rows
17c4762a1bSJed Brown
18c4762a1bSJed Brown Output parameter:
19c4762a1bSJed Brown sparsity - matrix sparsity pattern, typically computed using an ADOL-C function such as jac_pat
20c4762a1bSJed Brown */
PrintSparsity(MPI_Comm comm,PetscInt m,unsigned int ** sparsity)21d71ae5a4SJacob Faibussowitsch PetscErrorCode PrintSparsity(MPI_Comm comm, PetscInt m, unsigned int **sparsity)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown PetscFunctionBegin;
249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Sparsity pattern:\n"));
255f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < m; i++) {
269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n %2d: ", i));
2748a46eb9SPierre Jolivet for (PetscInt j = 1; j <= (PetscInt)sparsity[i][0]; j++) PetscCall(PetscPrintf(comm, " %2d ", sparsity[i][j]));
2860d4fc61SSatish Balay }
299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n\n"));
30*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
31c4762a1bSJed Brown }
32c4762a1bSJed Brown
33c4762a1bSJed Brown /*
34c4762a1bSJed Brown Generate a seed matrix defining the partition of columns of a matrix by a particular coloring,
35c4762a1bSJed Brown used for matrix compression
36c4762a1bSJed Brown
37c4762a1bSJed Brown Input parameter:
38c4762a1bSJed Brown iscoloring - the index set coloring to be used
39c4762a1bSJed Brown
40c4762a1bSJed Brown Output parameter:
41c4762a1bSJed Brown S - the resulting seed matrix
42c4762a1bSJed Brown
43c4762a1bSJed Brown Notes:
44c4762a1bSJed Brown Before calling GenerateSeedMatrix, Seed should be allocated as a logically 2d array with number of
45c4762a1bSJed Brown rows equal to the matrix to be compressed and number of columns equal to the number of colors used
46c4762a1bSJed Brown in iscoloring.
47c4762a1bSJed Brown */
GenerateSeedMatrix(ISColoring iscoloring,PetscScalar ** S)48d71ae5a4SJacob Faibussowitsch PetscErrorCode GenerateSeedMatrix(ISColoring iscoloring, PetscScalar **S)
49d71ae5a4SJacob Faibussowitsch {
50c4762a1bSJed Brown IS *is;
515f80ce2aSJacob Faibussowitsch PetscInt p, size;
52c4762a1bSJed Brown const PetscInt *indices;
53c4762a1bSJed Brown
54c4762a1bSJed Brown PetscFunctionBegin;
559566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &p, &is));
565f80ce2aSJacob Faibussowitsch for (PetscInt colour = 0; colour < p; colour++) {
579566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is[colour], &size));
589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[colour], &indices));
595f80ce2aSJacob Faibussowitsch for (PetscInt j = 0; j < size; j++) S[indices[j]][colour] = 1.;
609566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[colour], &indices));
61c4762a1bSJed Brown }
629566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(iscoloring, PETSC_USE_POINTER, &is));
63*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown
66c4762a1bSJed Brown /*
67c4762a1bSJed Brown Establish a look-up vector whose entries contain the colour used for that diagonal entry. Clearly
68c4762a1bSJed Brown we require the number of dependent and independent variables to be equal in this case.
69c4762a1bSJed Brown Input parameters:
70c4762a1bSJed Brown S - the seed matrix defining the coloring
71c4762a1bSJed Brown sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
72c4762a1bSJed Brown function, such as jac_pat or hess_pat
73c4762a1bSJed Brown m - the number of rows of Seed (and the matrix to be recovered)
74c4762a1bSJed Brown p - the number of colors used (also the number of columns in Seed)
75c4762a1bSJed Brown Output parameter:
76c4762a1bSJed Brown R - the recovery vector to be used for de-compression
77c4762a1bSJed Brown */
GenerateSeedMatrixPlusRecovery(ISColoring iscoloring,PetscScalar ** S,PetscScalar * R)78d71ae5a4SJacob Faibussowitsch PetscErrorCode GenerateSeedMatrixPlusRecovery(ISColoring iscoloring, PetscScalar **S, PetscScalar *R)
79d71ae5a4SJacob Faibussowitsch {
80c4762a1bSJed Brown IS *is;
81c4762a1bSJed Brown PetscInt p, size, colour, j;
82c4762a1bSJed Brown const PetscInt *indices;
83c4762a1bSJed Brown
84c4762a1bSJed Brown PetscFunctionBegin;
859566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &p, &is));
86c4762a1bSJed Brown for (colour = 0; colour < p; colour++) {
879566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is[colour], &size));
889566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[colour], &indices));
89c4762a1bSJed Brown for (j = 0; j < size; j++) {
90c4762a1bSJed Brown S[indices[j]][colour] = 1.;
91c4762a1bSJed Brown R[indices[j]] = colour;
92c4762a1bSJed Brown }
939566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is[colour], &indices));
94c4762a1bSJed Brown }
959566063dSJacob Faibussowitsch PetscCall(ISColoringRestoreIS(iscoloring, PETSC_USE_POINTER, &is));
96*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown Establish a look-up matrix whose entries contain the column coordinates of the corresponding entry
101c4762a1bSJed Brown in a matrix which has been compressed using the coloring defined by some seed matrix
102c4762a1bSJed Brown
103c4762a1bSJed Brown Input parameters:
104c4762a1bSJed Brown S - the seed matrix defining the coloring
105c4762a1bSJed Brown sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
106c4762a1bSJed Brown function, such as jac_pat or hess_pat
107c4762a1bSJed Brown m - the number of rows of Seed (and the matrix to be recovered)
108c4762a1bSJed Brown p - the number of colors used (also the number of columns in Seed)
109c4762a1bSJed Brown
110c4762a1bSJed Brown Output parameter:
111c4762a1bSJed Brown R - the recovery matrix to be used for de-compression
112c4762a1bSJed Brown */
GetRecoveryMatrix(PetscScalar ** S,unsigned int ** sparsity,PetscInt m,PetscInt p,PetscScalar ** R)113d71ae5a4SJacob Faibussowitsch PetscErrorCode GetRecoveryMatrix(PetscScalar **S, unsigned int **sparsity, PetscInt m, PetscInt p, PetscScalar **R)
114d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown PetscInt i, j, k, colour;
116c4762a1bSJed Brown
117c4762a1bSJed Brown PetscFunctionBegin;
118c4762a1bSJed Brown for (i = 0; i < m; i++) {
119c4762a1bSJed Brown for (colour = 0; colour < p; colour++) {
120c4762a1bSJed Brown R[i][colour] = -1.;
121c4762a1bSJed Brown for (k = 1; k <= (PetscInt)sparsity[i][0]; k++) {
122c4762a1bSJed Brown j = (PetscInt)sparsity[i][k];
123c4762a1bSJed Brown if (S[j][colour] == 1.) {
124c4762a1bSJed Brown R[i][colour] = j;
125c4762a1bSJed Brown break;
126c4762a1bSJed Brown }
127c4762a1bSJed Brown }
128c4762a1bSJed Brown }
129c4762a1bSJed Brown }
130*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown
133c4762a1bSJed Brown /*
134c4762a1bSJed Brown Recover the values of a sparse matrix from a compressed format and insert these into a matrix
135c4762a1bSJed Brown
136c4762a1bSJed Brown Input parameters:
137c4762a1bSJed Brown mode - use INSERT_VALUES or ADD_VALUES, as required
138c4762a1bSJed Brown m - number of rows of matrix.
139c4762a1bSJed Brown p - number of colors used in the compression of J (also the number of columns of R)
140c4762a1bSJed Brown R - recovery matrix to use in the decompression procedure
141c4762a1bSJed Brown C - compressed matrix to recover values from
142c4762a1bSJed Brown a - shift value for implicit problems (select NULL or unity for explicit problems)
143c4762a1bSJed Brown
144c4762a1bSJed Brown Output parameter:
145c4762a1bSJed Brown A - Mat to be populated with values from compressed matrix
146c4762a1bSJed Brown */
RecoverJacobian(Mat A,InsertMode mode,PetscInt m,PetscInt p,PetscScalar ** R,PetscScalar ** C,PetscReal * a)147d71ae5a4SJacob Faibussowitsch PetscErrorCode RecoverJacobian(Mat A, InsertMode mode, PetscInt m, PetscInt p, PetscScalar **R, PetscScalar **C, PetscReal *a)
148d71ae5a4SJacob Faibussowitsch {
149c4762a1bSJed Brown PetscFunctionBegin;
1505f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < m; i++) {
1515f80ce2aSJacob Faibussowitsch for (PetscInt colour = 0; colour < p; colour++) {
1525f80ce2aSJacob Faibussowitsch PetscInt j = (PetscInt)R[i][colour];
153c4762a1bSJed Brown if (j != -1) {
1545f80ce2aSJacob Faibussowitsch if (a) C[i][colour] *= *a;
1559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &C[i][colour], mode));
156c4762a1bSJed Brown }
157c4762a1bSJed Brown }
158c4762a1bSJed Brown }
159*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
160c4762a1bSJed Brown }
161c4762a1bSJed Brown
162c4762a1bSJed Brown /*
163c4762a1bSJed Brown Recover the values of the local portion of a sparse matrix from a compressed format and insert
164c4762a1bSJed Brown these into the local portion of a matrix
165c4762a1bSJed Brown
166c4762a1bSJed Brown Input parameters:
167c4762a1bSJed Brown mode - use INSERT_VALUES or ADD_VALUES, as required
168c4762a1bSJed Brown m - number of rows of matrix.
169c4762a1bSJed Brown p - number of colors used in the compression of J (also the number of columns of R)
170c4762a1bSJed Brown R - recovery matrix to use in the decompression procedure
171c4762a1bSJed Brown C - compressed matrix to recover values from
172c4762a1bSJed Brown a - shift value for implicit problems (select NULL or unity for explicit problems)
173c4762a1bSJed Brown
174c4762a1bSJed Brown Output parameter:
175c4762a1bSJed Brown A - Mat to be populated with values from compressed matrix
176c4762a1bSJed Brown */
RecoverJacobianLocal(Mat A,InsertMode mode,PetscInt m,PetscInt p,PetscScalar ** R,PetscScalar ** C,PetscReal * a)177d71ae5a4SJacob Faibussowitsch PetscErrorCode RecoverJacobianLocal(Mat A, InsertMode mode, PetscInt m, PetscInt p, PetscScalar **R, PetscScalar **C, PetscReal *a)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown PetscFunctionBegin;
1805f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < m; i++) {
1815f80ce2aSJacob Faibussowitsch for (PetscInt colour = 0; colour < p; colour++) {
1825f80ce2aSJacob Faibussowitsch PetscInt j = (PetscInt)R[i][colour];
183c4762a1bSJed Brown if (j != -1) {
1845f80ce2aSJacob Faibussowitsch if (a) C[i][colour] *= *a;
1859566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &C[i][colour], mode));
186c4762a1bSJed Brown }
187c4762a1bSJed Brown }
188c4762a1bSJed Brown }
189*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
192c4762a1bSJed Brown /*
193c4762a1bSJed Brown Recover the diagonal of the Jacobian from its compressed matrix format
194c4762a1bSJed Brown
195c4762a1bSJed Brown Input parameters:
196c4762a1bSJed Brown mode - use INSERT_VALUES or ADD_VALUES, as required
197c4762a1bSJed Brown m - number of rows of matrix.
198c4762a1bSJed Brown R - recovery vector to use in the decompression procedure
199c4762a1bSJed Brown C - compressed matrix to recover values from
200c4762a1bSJed Brown a - shift value for implicit problems (select NULL or unity for explicit problems)
201c4762a1bSJed Brown
202c4762a1bSJed Brown Output parameter:
203c4762a1bSJed Brown diag - Vec to be populated with values from compressed matrix
204c4762a1bSJed Brown */
RecoverDiagonal(Vec diag,InsertMode mode,PetscInt m,PetscScalar * R,PetscScalar ** C,PetscReal * a)205d71ae5a4SJacob Faibussowitsch PetscErrorCode RecoverDiagonal(Vec diag, InsertMode mode, PetscInt m, PetscScalar *R, PetscScalar **C, PetscReal *a)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown PetscFunctionBegin;
2085f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < m; i++) {
2095f80ce2aSJacob Faibussowitsch PetscInt colour = (PetscInt)R[i];
2105f80ce2aSJacob Faibussowitsch if (a) C[i][colour] *= *a;
2119566063dSJacob Faibussowitsch PetscCall(VecSetValues(diag, 1, &i, &C[i][colour], mode));
212c4762a1bSJed Brown }
213*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown
216c4762a1bSJed Brown /*
217c4762a1bSJed Brown Recover the local portion of the diagonal of the Jacobian from its compressed matrix format
218c4762a1bSJed Brown
219c4762a1bSJed Brown Input parameters:
220c4762a1bSJed Brown mode - use INSERT_VALUES or ADD_VALUES, as required
221c4762a1bSJed Brown m - number of rows of matrix.
222c4762a1bSJed Brown R - recovery vector to use in the decompression procedure
223c4762a1bSJed Brown C - compressed matrix to recover values from
224c4762a1bSJed Brown a - shift value for implicit problems (select NULL or unity for explicit problems)
225c4762a1bSJed Brown
226c4762a1bSJed Brown Output parameter:
227c4762a1bSJed Brown diag - Vec to be populated with values from compressed matrix
228c4762a1bSJed Brown */
RecoverDiagonalLocal(Vec diag,InsertMode mode,PetscInt m,PetscScalar * R,PetscScalar ** C,PetscReal * a)229d71ae5a4SJacob Faibussowitsch PetscErrorCode RecoverDiagonalLocal(Vec diag, InsertMode mode, PetscInt m, PetscScalar *R, PetscScalar **C, PetscReal *a)
230d71ae5a4SJacob Faibussowitsch {
231c4762a1bSJed Brown PetscFunctionBegin;
2325f80ce2aSJacob Faibussowitsch for (PetscInt i = 0; i < m; i++) {
2335f80ce2aSJacob Faibussowitsch PetscInt colour = (PetscInt)R[i];
2345f80ce2aSJacob Faibussowitsch if (a) C[i][colour] *= *a;
2359566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(diag, 1, &i, &C[i][colour], mode));
236c4762a1bSJed Brown }
237*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
238c4762a1bSJed Brown }
239