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