xref: /petsc/src/ts/tutorials/autodiff/adolc-utils/sparse.cxx (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
219371c9d4SSatish Balay PetscErrorCode PrintSparsity(MPI_Comm comm, PetscInt m, unsigned int **sparsity) {
22c4762a1bSJed Brown   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Sparsity pattern:\n"));
245f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < m; i++) {
259566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "\n %2d: ", i));
26*48a46eb9SPierre Jolivet     for (PetscInt j = 1; j <= (PetscInt)sparsity[i][0]; j++) PetscCall(PetscPrintf(comm, " %2d ", sparsity[i][j]));
2760d4fc61SSatish Balay   }
289566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "\n\n"));
29c4762a1bSJed Brown   PetscFunctionReturn(0);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown   Generate a seed matrix defining the partition of columns of a matrix by a particular coloring,
34c4762a1bSJed Brown   used for matrix compression
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   Input parameter:
37c4762a1bSJed Brown   iscoloring - the index set coloring to be used
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   Output parameter:
40c4762a1bSJed Brown   S          - the resulting seed matrix
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   Notes:
43c4762a1bSJed Brown   Before calling GenerateSeedMatrix, Seed should be allocated as a logically 2d array with number of
44c4762a1bSJed Brown   rows equal to the matrix to be compressed and number of columns equal to the number of colors used
45c4762a1bSJed Brown   in iscoloring.
46c4762a1bSJed Brown */
479371c9d4SSatish Balay PetscErrorCode GenerateSeedMatrix(ISColoring iscoloring, PetscScalar **S) {
48c4762a1bSJed Brown   IS             *is;
495f80ce2aSJacob Faibussowitsch   PetscInt        p, size;
50c4762a1bSJed Brown   const PetscInt *indices;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &p, &is));
545f80ce2aSJacob Faibussowitsch   for (PetscInt colour = 0; colour < p; colour++) {
559566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is[colour], &size));
569566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[colour], &indices));
575f80ce2aSJacob Faibussowitsch     for (PetscInt j = 0; j < size; j++) S[indices[j]][colour] = 1.;
589566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[colour], &indices));
59c4762a1bSJed Brown   }
609566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_USE_POINTER, &is));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*
65c4762a1bSJed Brown   Establish a look-up vector whose entries contain the colour used for that diagonal entry. Clearly
66c4762a1bSJed Brown   we require the number of dependent and independent variables to be equal in this case.
67c4762a1bSJed Brown   Input parameters:
68c4762a1bSJed Brown   S        - the seed matrix defining the coloring
69c4762a1bSJed Brown   sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
70c4762a1bSJed Brown              function, such as jac_pat or hess_pat
71c4762a1bSJed Brown   m        - the number of rows of Seed (and the matrix to be recovered)
72c4762a1bSJed Brown   p        - the number of colors used (also the number of columns in Seed)
73c4762a1bSJed Brown   Output parameter:
74c4762a1bSJed Brown   R        - the recovery vector to be used for de-compression
75c4762a1bSJed Brown */
769371c9d4SSatish Balay PetscErrorCode GenerateSeedMatrixPlusRecovery(ISColoring iscoloring, PetscScalar **S, PetscScalar *R) {
77c4762a1bSJed Brown   IS             *is;
78c4762a1bSJed Brown   PetscInt        p, size, colour, j;
79c4762a1bSJed Brown   const PetscInt *indices;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &p, &is));
83c4762a1bSJed Brown   for (colour = 0; colour < p; colour++) {
849566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(is[colour], &size));
859566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is[colour], &indices));
86c4762a1bSJed Brown     for (j = 0; j < size; j++) {
87c4762a1bSJed Brown       S[indices[j]][colour] = 1.;
88c4762a1bSJed Brown       R[indices[j]]         = colour;
89c4762a1bSJed Brown     }
909566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is[colour], &indices));
91c4762a1bSJed Brown   }
929566063dSJacob Faibussowitsch   PetscCall(ISColoringRestoreIS(iscoloring, PETSC_USE_POINTER, &is));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /*
97c4762a1bSJed Brown   Establish a look-up matrix whose entries contain the column coordinates of the corresponding entry
98c4762a1bSJed Brown   in a matrix which has been compressed using the coloring defined by some seed matrix
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   Input parameters:
101c4762a1bSJed Brown   S        - the seed matrix defining the coloring
102c4762a1bSJed Brown   sparsity - the sparsity pattern of the matrix to be recovered, typically computed using an ADOL-C
103c4762a1bSJed Brown              function, such as jac_pat or hess_pat
104c4762a1bSJed Brown   m        - the number of rows of Seed (and the matrix to be recovered)
105c4762a1bSJed Brown   p        - the number of colors used (also the number of columns in Seed)
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   Output parameter:
108c4762a1bSJed Brown   R        - the recovery matrix to be used for de-compression
109c4762a1bSJed Brown */
1109371c9d4SSatish Balay PetscErrorCode GetRecoveryMatrix(PetscScalar **S, unsigned int **sparsity, PetscInt m, PetscInt p, PetscScalar **R) {
111c4762a1bSJed Brown   PetscInt i, j, k, colour;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBegin;
114c4762a1bSJed Brown   for (i = 0; i < m; i++) {
115c4762a1bSJed Brown     for (colour = 0; colour < p; colour++) {
116c4762a1bSJed Brown       R[i][colour] = -1.;
117c4762a1bSJed Brown       for (k = 1; k <= (PetscInt)sparsity[i][0]; k++) {
118c4762a1bSJed Brown         j = (PetscInt)sparsity[i][k];
119c4762a1bSJed Brown         if (S[j][colour] == 1.) {
120c4762a1bSJed Brown           R[i][colour] = j;
121c4762a1bSJed Brown           break;
122c4762a1bSJed Brown         }
123c4762a1bSJed Brown       }
124c4762a1bSJed Brown     }
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown   PetscFunctionReturn(0);
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown /*
130c4762a1bSJed Brown   Recover the values of a sparse matrix from a compressed format and insert these into a matrix
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   Input parameters:
133c4762a1bSJed Brown   mode - use INSERT_VALUES or ADD_VALUES, as required
134c4762a1bSJed Brown   m    - number of rows of matrix.
135c4762a1bSJed Brown   p    - number of colors used in the compression of J (also the number of columns of R)
136c4762a1bSJed Brown   R    - recovery matrix to use in the decompression procedure
137c4762a1bSJed Brown   C    - compressed matrix to recover values from
138c4762a1bSJed Brown   a    - shift value for implicit problems (select NULL or unity for explicit problems)
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   Output parameter:
141c4762a1bSJed Brown   A    - Mat to be populated with values from compressed matrix
142c4762a1bSJed Brown */
1439371c9d4SSatish Balay PetscErrorCode RecoverJacobian(Mat A, InsertMode mode, PetscInt m, PetscInt p, PetscScalar **R, PetscScalar **C, PetscReal *a) {
144c4762a1bSJed Brown   PetscFunctionBegin;
1455f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < m; i++) {
1465f80ce2aSJacob Faibussowitsch     for (PetscInt colour = 0; colour < p; colour++) {
1475f80ce2aSJacob Faibussowitsch       PetscInt j = (PetscInt)R[i][colour];
148c4762a1bSJed Brown       if (j != -1) {
1495f80ce2aSJacob Faibussowitsch         if (a) C[i][colour] *= *a;
1509566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &i, 1, &j, &C[i][colour], mode));
151c4762a1bSJed Brown       }
152c4762a1bSJed Brown     }
153c4762a1bSJed Brown   }
154c4762a1bSJed Brown   PetscFunctionReturn(0);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown /*
158c4762a1bSJed Brown   Recover the values of the local portion of a sparse matrix from a compressed format and insert
159c4762a1bSJed Brown   these into the local portion of a matrix
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   Input parameters:
162c4762a1bSJed Brown   mode - use INSERT_VALUES or ADD_VALUES, as required
163c4762a1bSJed Brown   m    - number of rows of matrix.
164c4762a1bSJed Brown   p    - number of colors used in the compression of J (also the number of columns of R)
165c4762a1bSJed Brown   R    - recovery matrix to use in the decompression procedure
166c4762a1bSJed Brown   C    - compressed matrix to recover values from
167c4762a1bSJed Brown   a    - shift value for implicit problems (select NULL or unity for explicit problems)
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   Output parameter:
170c4762a1bSJed Brown   A    - Mat to be populated with values from compressed matrix
171c4762a1bSJed Brown */
1729371c9d4SSatish Balay PetscErrorCode RecoverJacobianLocal(Mat A, InsertMode mode, PetscInt m, PetscInt p, PetscScalar **R, PetscScalar **C, PetscReal *a) {
173c4762a1bSJed Brown   PetscFunctionBegin;
1745f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < m; i++) {
1755f80ce2aSJacob Faibussowitsch     for (PetscInt colour = 0; colour < p; colour++) {
1765f80ce2aSJacob Faibussowitsch       PetscInt j = (PetscInt)R[i][colour];
177c4762a1bSJed Brown       if (j != -1) {
1785f80ce2aSJacob Faibussowitsch         if (a) C[i][colour] *= *a;
1799566063dSJacob Faibussowitsch         PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &C[i][colour], mode));
180c4762a1bSJed Brown       }
181c4762a1bSJed Brown     }
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown   PetscFunctionReturn(0);
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /*
187c4762a1bSJed Brown   Recover the diagonal of the Jacobian from its compressed matrix format
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   Input parameters:
190c4762a1bSJed Brown   mode - use INSERT_VALUES or ADD_VALUES, as required
191c4762a1bSJed Brown   m    - number of rows of matrix.
192c4762a1bSJed Brown   R    - recovery vector to use in the decompression procedure
193c4762a1bSJed Brown   C    - compressed matrix to recover values from
194c4762a1bSJed Brown   a    - shift value for implicit problems (select NULL or unity for explicit problems)
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   Output parameter:
197c4762a1bSJed Brown   diag - Vec to be populated with values from compressed matrix
198c4762a1bSJed Brown */
1999371c9d4SSatish Balay PetscErrorCode RecoverDiagonal(Vec diag, InsertMode mode, PetscInt m, PetscScalar *R, PetscScalar **C, PetscReal *a) {
200c4762a1bSJed Brown   PetscFunctionBegin;
2015f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < m; i++) {
2025f80ce2aSJacob Faibussowitsch     PetscInt colour = (PetscInt)R[i];
2035f80ce2aSJacob Faibussowitsch     if (a) C[i][colour] *= *a;
2049566063dSJacob Faibussowitsch     PetscCall(VecSetValues(diag, 1, &i, &C[i][colour], mode));
205c4762a1bSJed Brown   }
206c4762a1bSJed Brown   PetscFunctionReturn(0);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /*
210c4762a1bSJed Brown   Recover the local portion of the diagonal of the Jacobian from its compressed matrix format
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   Input parameters:
213c4762a1bSJed Brown   mode - use INSERT_VALUES or ADD_VALUES, as required
214c4762a1bSJed Brown   m    - number of rows of matrix.
215c4762a1bSJed Brown   R    - recovery vector to use in the decompression procedure
216c4762a1bSJed Brown   C    - compressed matrix to recover values from
217c4762a1bSJed Brown   a    - shift value for implicit problems (select NULL or unity for explicit problems)
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   Output parameter:
220c4762a1bSJed Brown   diag - Vec to be populated with values from compressed matrix
221c4762a1bSJed Brown */
2229371c9d4SSatish Balay PetscErrorCode RecoverDiagonalLocal(Vec diag, InsertMode mode, PetscInt m, PetscScalar *R, PetscScalar **C, PetscReal *a) {
223c4762a1bSJed Brown   PetscFunctionBegin;
2245f80ce2aSJacob Faibussowitsch   for (PetscInt i = 0; i < m; i++) {
2255f80ce2aSJacob Faibussowitsch     PetscInt colour = (PetscInt)R[i];
2265f80ce2aSJacob Faibussowitsch     if (a) C[i][colour] *= *a;
2279566063dSJacob Faibussowitsch     PetscCall(VecSetValuesLocal(diag, 1, &i, &C[i][colour], mode));
228c4762a1bSJed Brown   }
229c4762a1bSJed Brown   PetscFunctionReturn(0);
230c4762a1bSJed Brown }
231