1c4762a1bSJed Brown #include "contexts.cxx" 2c4762a1bSJed Brown #include "sparse.cxx" 3c4762a1bSJed Brown #include "init.cxx" 4c4762a1bSJed Brown #include <adolc/drivers/drivers.h> 5c4762a1bSJed Brown #include <adolc/interfaces.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 9c4762a1bSJed Brown 10c4762a1bSJed Brown For documentation on ADOL-C, see 11c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown /* -------------------------------------------------------------------------------- 15c4762a1bSJed Brown Drivers for RHSJacobian and IJacobian 16c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown Compute Jacobian for explicit TS in compressed format and recover from this, using 20c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 21c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 22c4762a1bSJed Brown 23c4762a1bSJed Brown Input parameters: 24c4762a1bSJed Brown tag - tape identifier 25c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 26c4762a1bSJed Brown ctx - ADOL-C context, as defined above 27c4762a1bSJed Brown 28c4762a1bSJed Brown Output parameter: 29c4762a1bSJed Brown A - Mat object corresponding to Jacobian 30c4762a1bSJed Brown */ 319371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeRHSJacobian(PetscInt tag, Mat A, const PetscScalar *u_vec, void *ctx) { 32c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 33c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 34c4762a1bSJed Brown PetscScalar **J; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBegin; 379566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 389371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J); 399371c9d4SSatish Balay else jacobian(tag, m, n, u_vec, J); 40c4762a1bSJed Brown if (adctx->sparse) { 419566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 42c4762a1bSJed Brown } else { 43c4762a1bSJed Brown for (i = 0; i < m; i++) { 44c4762a1bSJed Brown for (j = 0; j < n; j++) { 45*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown } 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown Compute Jacobian for explicit TS in compressed format and recover from this, using 57c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 58c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 59c4762a1bSJed Brown 60c4762a1bSJed Brown Input parameters: 61c4762a1bSJed Brown tag - tape identifier 62c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 63c4762a1bSJed Brown ctx - ADOL-C context, as defined above 64c4762a1bSJed Brown 65c4762a1bSJed Brown Output parameter: 66c4762a1bSJed Brown A - Mat object corresponding to Jacobian 67c4762a1bSJed Brown */ 689371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeRHSJacobianLocal(PetscInt tag, Mat A, const PetscScalar *u_vec, void *ctx) { 69c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 70c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 71c4762a1bSJed Brown PetscScalar **J; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBegin; 749566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 759371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J); 769371c9d4SSatish Balay else jacobian(tag, m, n, u_vec, J); 77c4762a1bSJed Brown if (adctx->sparse) { 789566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 79c4762a1bSJed Brown } else { 80c4762a1bSJed Brown for (i = 0; i < m; i++) { 81c4762a1bSJed Brown for (j = 0; j < n; j++) { 82*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown } 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown PetscFunctionReturn(0); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown Compute Jacobian for implicit TS in compressed format and recover from this, using 94c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 95c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 96c4762a1bSJed Brown 97c4762a1bSJed Brown Input parameters: 98c4762a1bSJed Brown tag1 - tape identifier for df/dx part 99c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 100c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 101c4762a1bSJed Brown ctx - ADOL-C context, as defined above 102c4762a1bSJed Brown 103c4762a1bSJed Brown Output parameter: 104c4762a1bSJed Brown A - Mat object corresponding to Jacobian 105c4762a1bSJed Brown */ 1069371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeIJacobian(PetscInt tag1, PetscInt tag2, Mat A, const PetscScalar *u_vec, PetscReal a, void *ctx) { 107c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 108c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 109c4762a1bSJed Brown PetscScalar **J; 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* dF/dx part */ 1159371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J); 1169371c9d4SSatish Balay else jacobian(tag1, m, n, u_vec, J); 1179566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 118c4762a1bSJed Brown if (adctx->sparse) { 1199566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 120c4762a1bSJed Brown } else { 121c4762a1bSJed Brown for (i = 0; i < m; i++) { 122c4762a1bSJed Brown for (j = 0; j < n; j++) { 123*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126c4762a1bSJed Brown } 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* a * dF/d(xdot) part */ 1319371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J); 1329371c9d4SSatish Balay else jacobian(tag2, m, n, u_vec, J); 133c4762a1bSJed Brown if (adctx->sparse) { 1349566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A, ADD_VALUES, m, p, adctx->Rec, J, &a)); 135c4762a1bSJed Brown } else { 136c4762a1bSJed Brown for (i = 0; i < m; i++) { 137c4762a1bSJed Brown for (j = 0; j < n; j++) { 138c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 139c4762a1bSJed Brown J[i][j] *= a; 1409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], ADD_VALUES)); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1479566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 148c4762a1bSJed Brown PetscFunctionReturn(0); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* 152c4762a1bSJed Brown Compute Jacobian for implicit TS in the special case where it is 153c4762a1bSJed Brown known that the mass matrix is simply the identity. i.e. We have 154c4762a1bSJed Brown a problem of the form 155c4762a1bSJed Brown du/dt = F(u). 156c4762a1bSJed Brown 157c4762a1bSJed Brown Input parameters: 158c4762a1bSJed Brown tag - tape identifier for df/dx part 159c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 160c4762a1bSJed Brown ctx - ADOL-C context, as defined above 161c4762a1bSJed Brown 162c4762a1bSJed Brown Output parameter: 163c4762a1bSJed Brown A - Mat object corresponding to Jacobian 164c4762a1bSJed Brown */ 1659371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeIJacobianIDMass(PetscInt tag, Mat A, PetscScalar *u_vec, PetscReal a, void *ctx) { 166c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 167c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 168c4762a1bSJed Brown PetscScalar **J; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBegin; 1719566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* dF/dx part */ 1749371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J); 1759371c9d4SSatish Balay else jacobian(tag, m, n, u_vec, J); 1769566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 177c4762a1bSJed Brown if (adctx->sparse) { 1789566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 179c4762a1bSJed Brown } else { 180c4762a1bSJed Brown for (i = 0; i < m; i++) { 181c4762a1bSJed Brown for (j = 0; j < n; j++) { 182*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 1869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1889566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* a * dF/d(xdot) part */ 1919566063dSJacob Faibussowitsch PetscCall(MatShift(A, a)); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown Compute local portion of Jacobian for implicit TS in compressed format and recover from this, using 197c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 198c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 199c4762a1bSJed Brown 200c4762a1bSJed Brown Input parameters: 201c4762a1bSJed Brown tag1 - tape identifier for df/dx part 202c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 203c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 204c4762a1bSJed Brown ctx - ADOL-C context, as defined above 205c4762a1bSJed Brown 206c4762a1bSJed Brown Output parameter: 207c4762a1bSJed Brown A - Mat object corresponding to Jacobian 208c4762a1bSJed Brown */ 2099371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeIJacobianLocal(PetscInt tag1, PetscInt tag2, Mat A, PetscScalar *u_vec, PetscReal a, void *ctx) { 210c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 211c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 212c4762a1bSJed Brown PetscScalar **J; 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscFunctionBegin; 2159566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* dF/dx part */ 2189371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J); 2199371c9d4SSatish Balay else jacobian(tag1, m, n, u_vec, J); 220c4762a1bSJed Brown if (adctx->sparse) { 2219566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 222c4762a1bSJed Brown } else { 223c4762a1bSJed Brown for (i = 0; i < m; i++) { 224c4762a1bSJed Brown for (j = 0; j < n; j++) { 225*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 2299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* a * dF/d(xdot) part */ 2339371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J); 2349371c9d4SSatish Balay else jacobian(tag2, m, n, u_vec, J); 235c4762a1bSJed Brown if (adctx->sparse) { 2369566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A, ADD_VALUES, m, p, adctx->Rec, J, &a)); 237c4762a1bSJed Brown } else { 238c4762a1bSJed Brown for (i = 0; i < m; i++) { 239c4762a1bSJed Brown for (j = 0; j < n; j++) { 240c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 241c4762a1bSJed Brown J[i][j] *= a; 2429566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], ADD_VALUES)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 2479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2499566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 250c4762a1bSJed Brown PetscFunctionReturn(0); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown Compute local portion of Jacobian for implicit TS in the special case where it is 255c4762a1bSJed Brown known that the mass matrix is simply the identity. i.e. We have 256c4762a1bSJed Brown a problem of the form 257c4762a1bSJed Brown du/dt = F(u). 258c4762a1bSJed Brown 259c4762a1bSJed Brown Input parameters: 260c4762a1bSJed Brown tag - tape identifier for df/dx part 261c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 262c4762a1bSJed Brown ctx - ADOL-C context, as defined above 263c4762a1bSJed Brown 264c4762a1bSJed Brown Output parameter: 265c4762a1bSJed Brown A - Mat object corresponding to Jacobian 266c4762a1bSJed Brown */ 2679371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeIJacobianLocalIDMass(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscReal a, void *ctx) { 268c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 269c4762a1bSJed Brown PetscInt i, j, m = adctx->m, n = adctx->n, p = adctx->p; 270c4762a1bSJed Brown PetscScalar **J; 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* dF/dx part */ 2769371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J); 2779371c9d4SSatish Balay else jacobian(tag, m, n, u_vec, J); 278c4762a1bSJed Brown if (adctx->sparse) { 2799566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL)); 280c4762a1bSJed Brown } else { 281c4762a1bSJed Brown for (i = 0; i < m; i++) { 282c4762a1bSJed Brown for (j = 0; j < n; j++) { 283*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 2879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2899566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* a * dF/d(xdot) part */ 2929566063dSJacob Faibussowitsch PetscCall(MatShift(A, a)); 293c4762a1bSJed Brown PetscFunctionReturn(0); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* -------------------------------------------------------------------------------- 297c4762a1bSJed Brown Drivers for Jacobian w.r.t. a parameter 298c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* 301c4762a1bSJed Brown Compute Jacobian w.r.t a parameter for explicit TS. 302c4762a1bSJed Brown 303c4762a1bSJed Brown Input parameters: 304c4762a1bSJed Brown tag - tape identifier 305c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 306c4762a1bSJed Brown params - the parameters w.r.t. which we differentiate 307c4762a1bSJed Brown ctx - ADOL-C context, as defined above 308c4762a1bSJed Brown 309c4762a1bSJed Brown Output parameter: 310c4762a1bSJed Brown A - Mat object corresponding to Jacobian 311c4762a1bSJed Brown */ 3129371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeRHSJacobianP(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscScalar *params, void *ctx) { 313c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 314c4762a1bSJed Brown PetscInt i, j = 0, m = adctx->m, n = adctx->n, p = adctx->num_params; 315c4762a1bSJed Brown PetscScalar **J, *concat, **S; 316c4762a1bSJed Brown 317c4762a1bSJed Brown PetscFunctionBegin; 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* Allocate memory and concatenate independent variable values with parameter */ 3209566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 3219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + p, &concat)); 3229566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(n + p, p, &S)); 3239566063dSJacob Faibussowitsch PetscCall(Subidentity(p, n, S)); 324c4762a1bSJed Brown for (i = 0; i < n; i++) concat[i] = u_vec[i]; 325c4762a1bSJed Brown for (i = 0; i < p; i++) concat[n + i] = params[i]; 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* Propagate the appropriate seed matrix through the forward mode of AD */ 328c4762a1bSJed Brown fov_forward(tag, m, n + p, p, concat, S, NULL, J); 3299566063dSJacob Faibussowitsch PetscCall(AdolcFree2(S)); 3309566063dSJacob Faibussowitsch PetscCall(PetscFree(concat)); 331c4762a1bSJed Brown 332c4762a1bSJed Brown /* Set matrix values */ 333c4762a1bSJed Brown for (i = 0; i < m; i++) { 334c4762a1bSJed Brown for (j = 0; j < p; j++) { 335*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } 3389566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 3399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 341c4762a1bSJed Brown PetscFunctionReturn(0); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* 345c4762a1bSJed Brown Compute local portion of Jacobian w.r.t a parameter for explicit TS. 346c4762a1bSJed Brown 347c4762a1bSJed Brown Input parameters: 348c4762a1bSJed Brown tag - tape identifier 349c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 350c4762a1bSJed Brown params - the parameters w.r.t. which we differentiate 351c4762a1bSJed Brown ctx - ADOL-C context, as defined above 352c4762a1bSJed Brown 353c4762a1bSJed Brown Output parameter: 354c4762a1bSJed Brown A - Mat object corresponding to Jacobian 355c4762a1bSJed Brown */ 3569371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeRHSJacobianPLocal(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscScalar *params, void *ctx) { 357c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 358c4762a1bSJed Brown PetscInt i, j = 0, m = adctx->m, n = adctx->n, p = adctx->num_params; 359c4762a1bSJed Brown PetscScalar **J, *concat, **S; 360c4762a1bSJed Brown 361c4762a1bSJed Brown PetscFunctionBegin; 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* Allocate memory and concatenate independent variable values with parameter */ 3649566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 3659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + p, &concat)); 3669566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(n + p, p, &S)); 3679566063dSJacob Faibussowitsch PetscCall(Subidentity(p, n, S)); 368c4762a1bSJed Brown for (i = 0; i < n; i++) concat[i] = u_vec[i]; 369c4762a1bSJed Brown for (i = 0; i < p; i++) concat[n + i] = params[i]; 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* Propagate the appropriate seed matrix through the forward mode of AD */ 372c4762a1bSJed Brown fov_forward(tag, m, n + p, p, concat, S, NULL, J); 3739566063dSJacob Faibussowitsch PetscCall(AdolcFree2(S)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(concat)); 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* Set matrix values */ 377c4762a1bSJed Brown for (i = 0; i < m; i++) { 378c4762a1bSJed Brown for (j = 0; j < p; j++) { 379*48a46eb9SPierre Jolivet if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } 3829566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 3839566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3849566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 385c4762a1bSJed Brown PetscFunctionReturn(0); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* -------------------------------------------------------------------------------- 389c4762a1bSJed Brown Drivers for Jacobian diagonal 390c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 391c4762a1bSJed Brown 392c4762a1bSJed Brown /* 393c4762a1bSJed Brown Compute local portion of Jacobian diagonal for implicit TS in compressed format and recover 394c4762a1bSJed Brown from this, using precomputed seed matrix and recovery vector. 395c4762a1bSJed Brown 396c4762a1bSJed Brown Input parameters: 397c4762a1bSJed Brown tag1 - tape identifier for df/dx part 398c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 399c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 400c4762a1bSJed Brown ctx - ADOL-C context, as defined above 401c4762a1bSJed Brown 402c4762a1bSJed Brown Output parameter: 403c4762a1bSJed Brown diag - Vec object corresponding to Jacobian diagonal 404c4762a1bSJed Brown */ 4059371c9d4SSatish Balay PetscErrorCode PetscAdolcComputeIJacobianAndDiagonalLocal(PetscInt tag1, PetscInt tag2, Vec diag, PetscScalar *u_vec, PetscReal a, void *ctx) { 406c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx *)ctx; 407c4762a1bSJed Brown PetscInt i, m = adctx->m, n = adctx->n, p = adctx->p; 408c4762a1bSJed Brown PetscScalar **J; 409c4762a1bSJed Brown 410c4762a1bSJed Brown PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m, p, &J)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* dF/dx part */ 4149371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J); 4159371c9d4SSatish Balay else jacobian(tag1, m, n, u_vec, J); 416c4762a1bSJed Brown if (adctx->sparse) { 4179566063dSJacob Faibussowitsch PetscCall(RecoverDiagonalLocal(diag, INSERT_VALUES, m, adctx->rec, J, NULL)); 418c4762a1bSJed Brown } else { 419c4762a1bSJed Brown for (i = 0; i < m; i++) { 420*48a46eb9SPierre Jolivet if (fabs(J[i][i]) > 1.e-16) PetscCall(VecSetValuesLocal(diag, 1, &i, &J[i][i], INSERT_VALUES)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown } 4239566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(diag)); 4249566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(diag)); 425c4762a1bSJed Brown 426c4762a1bSJed Brown /* a * dF/d(xdot) part */ 4279371c9d4SSatish Balay if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J); 4289371c9d4SSatish Balay else jacobian(tag2, m, n, u_vec, J); 429c4762a1bSJed Brown if (adctx->sparse) { 4309566063dSJacob Faibussowitsch PetscCall(RecoverDiagonalLocal(diag, ADD_VALUES, m, adctx->rec, J, NULL)); 431c4762a1bSJed Brown } else { 432c4762a1bSJed Brown for (i = 0; i < m; i++) { 433c4762a1bSJed Brown if (fabs(J[i][i]) > 1.e-16) { 434c4762a1bSJed Brown J[i][i] *= a; 4359566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(diag, 1, &i, &J[i][i], ADD_VALUES)); 436c4762a1bSJed Brown } 437c4762a1bSJed Brown } 438c4762a1bSJed Brown } 4399566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(diag)); 4409566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(diag)); 4419566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 442c4762a1bSJed Brown PetscFunctionReturn(0); 443c4762a1bSJed Brown } 444