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