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