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 */ 31a8c08197SHong Zhang PetscErrorCode PetscAdolcComputeRHSJacobian(PetscInt tag,Mat A,const PetscScalar *u_vec,void *ctx) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 34c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 35c4762a1bSJed Brown PetscScalar **J; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBegin; 38*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 39c4762a1bSJed Brown if (adctx->Seed) 40c4762a1bSJed Brown fov_forward(tag,m,n,p,u_vec,adctx->Seed,NULL,J); 41c4762a1bSJed Brown else 42c4762a1bSJed Brown jacobian(tag,m,n,u_vec,J); 43c4762a1bSJed Brown if (adctx->sparse) { 44*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 45c4762a1bSJed Brown } else { 46c4762a1bSJed Brown for (i=0; i<m; i++) { 47c4762a1bSJed Brown for (j=0; j<n; j++) { 48c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 49*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 54*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 55*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 56*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 57c4762a1bSJed Brown PetscFunctionReturn(0); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* 61c4762a1bSJed Brown Compute Jacobian for explicit TS in compressed format and recover from this, using 62c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 63c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 64c4762a1bSJed Brown 65c4762a1bSJed Brown Input parameters: 66c4762a1bSJed Brown tag - tape identifier 67c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 68c4762a1bSJed Brown ctx - ADOL-C context, as defined above 69c4762a1bSJed Brown 70c4762a1bSJed Brown Output parameter: 71c4762a1bSJed Brown A - Mat object corresponding to Jacobian 72c4762a1bSJed Brown */ 73a8c08197SHong Zhang PetscErrorCode PetscAdolcComputeRHSJacobianLocal(PetscInt tag,Mat A,const PetscScalar *u_vec,void *ctx) 74c4762a1bSJed Brown { 75c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 76c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 77c4762a1bSJed Brown PetscScalar **J; 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscFunctionBegin; 80*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 81c4762a1bSJed Brown if (adctx->Seed) 82c4762a1bSJed Brown fov_forward(tag,m,n,p,u_vec,adctx->Seed,NULL,J); 83c4762a1bSJed Brown else 84c4762a1bSJed Brown jacobian(tag,m,n,u_vec,J); 85c4762a1bSJed Brown if (adctx->sparse) { 86*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 87c4762a1bSJed Brown } else { 88c4762a1bSJed Brown for (i=0; i<m; i++) { 89c4762a1bSJed Brown for (j=0; j<n; j++) { 90c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 91*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown } 94c4762a1bSJed Brown } 95c4762a1bSJed Brown } 96*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 97*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 98*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Compute Jacobian for implicit TS in compressed format and recover from this, using 104c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 105c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 106c4762a1bSJed Brown 107c4762a1bSJed Brown Input parameters: 108c4762a1bSJed Brown tag1 - tape identifier for df/dx part 109c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 110c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 111c4762a1bSJed Brown ctx - ADOL-C context, as defined above 112c4762a1bSJed Brown 113c4762a1bSJed Brown Output parameter: 114c4762a1bSJed Brown A - Mat object corresponding to Jacobian 115c4762a1bSJed Brown */ 116410585f6SHong Zhang PetscErrorCode PetscAdolcComputeIJacobian(PetscInt tag1,PetscInt tag2,Mat A,const PetscScalar *u_vec,PetscReal a,void *ctx) 117c4762a1bSJed Brown { 118c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 119c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 120c4762a1bSJed Brown PetscScalar **J; 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscFunctionBegin; 123*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* dF/dx part */ 126c4762a1bSJed Brown if (adctx->Seed) 127c4762a1bSJed Brown fov_forward(tag1,m,n,p,u_vec,adctx->Seed,NULL,J); 128c4762a1bSJed Brown else 129c4762a1bSJed Brown jacobian(tag1,m,n,u_vec,J); 130*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 131c4762a1bSJed Brown if (adctx->sparse) { 132*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 133c4762a1bSJed Brown } else { 134c4762a1bSJed Brown for (i=0; i<m; i++) { 135c4762a1bSJed Brown for (j=0; j<n; j++) { 136c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 137*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 143*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* a * dF/d(xdot) part */ 146c4762a1bSJed Brown if (adctx->Seed) 147c4762a1bSJed Brown fov_forward(tag2,m,n,p,u_vec,adctx->Seed,NULL,J); 148c4762a1bSJed Brown else 149c4762a1bSJed Brown jacobian(tag2,m,n,u_vec,J); 150c4762a1bSJed Brown if (adctx->sparse) { 151*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A,ADD_VALUES,m,p,adctx->Rec,J,&a)); 152c4762a1bSJed Brown } else { 153c4762a1bSJed Brown for (i=0; i<m; i++) { 154c4762a1bSJed Brown for (j=0; j<n; j++) { 155c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 156c4762a1bSJed Brown J[i][j] *= a; 157*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],ADD_VALUES)); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 163*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 164*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 165c4762a1bSJed Brown PetscFunctionReturn(0); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* 169c4762a1bSJed Brown Compute Jacobian for implicit TS in the special case where it is 170c4762a1bSJed Brown known that the mass matrix is simply the identity. i.e. We have 171c4762a1bSJed Brown a problem of the form 172c4762a1bSJed Brown du/dt = F(u). 173c4762a1bSJed Brown 174c4762a1bSJed Brown Input parameters: 175c4762a1bSJed Brown tag - tape identifier for df/dx part 176c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 177c4762a1bSJed Brown ctx - ADOL-C context, as defined above 178c4762a1bSJed Brown 179c4762a1bSJed Brown Output parameter: 180c4762a1bSJed Brown A - Mat object corresponding to Jacobian 181c4762a1bSJed Brown */ 182c4762a1bSJed Brown PetscErrorCode PetscAdolcComputeIJacobianIDMass(PetscInt tag,Mat A,PetscScalar *u_vec,PetscReal a,void *ctx) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 185c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 186c4762a1bSJed Brown PetscScalar **J; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBegin; 189*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* dF/dx part */ 192c4762a1bSJed Brown if (adctx->Seed) 193c4762a1bSJed Brown fov_forward(tag,m,n,p,u_vec,adctx->Seed,NULL,J); 194c4762a1bSJed Brown else 195c4762a1bSJed Brown jacobian(tag,m,n,u_vec,J); 196*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 197c4762a1bSJed Brown if (adctx->sparse) { 198*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobian(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 199c4762a1bSJed Brown } else { 200c4762a1bSJed Brown for (i=0; i<m; i++) { 201c4762a1bSJed Brown for (j=0; j<n; j++) { 202c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 203*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown } 206c4762a1bSJed Brown } 207c4762a1bSJed Brown } 208*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 209*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 210*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* a * dF/d(xdot) part */ 213*9566063dSJacob Faibussowitsch PetscCall(MatShift(A,a)); 214c4762a1bSJed Brown PetscFunctionReturn(0); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* 218c4762a1bSJed Brown Compute local portion of Jacobian for implicit TS in compressed format and recover from this, using 219c4762a1bSJed Brown precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is 220c4762a1bSJed Brown assembled (not recommended for non-toy problems!). 221c4762a1bSJed Brown 222c4762a1bSJed Brown Input parameters: 223c4762a1bSJed Brown tag1 - tape identifier for df/dx part 224c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 225c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 226c4762a1bSJed Brown ctx - ADOL-C context, as defined above 227c4762a1bSJed Brown 228c4762a1bSJed Brown Output parameter: 229c4762a1bSJed Brown A - Mat object corresponding to Jacobian 230c4762a1bSJed Brown */ 231c4762a1bSJed Brown PetscErrorCode PetscAdolcComputeIJacobianLocal(PetscInt tag1,PetscInt tag2,Mat A,PetscScalar *u_vec,PetscReal a,void *ctx) 232c4762a1bSJed Brown { 233c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 234c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 235c4762a1bSJed Brown PetscScalar **J; 236c4762a1bSJed Brown 237c4762a1bSJed Brown PetscFunctionBegin; 238*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* dF/dx part */ 241c4762a1bSJed Brown if (adctx->Seed) 242c4762a1bSJed Brown fov_forward(tag1,m,n,p,u_vec,adctx->Seed,NULL,J); 243c4762a1bSJed Brown else 244c4762a1bSJed Brown jacobian(tag1,m,n,u_vec,J); 245c4762a1bSJed Brown if (adctx->sparse) { 246*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 247c4762a1bSJed Brown } else { 248c4762a1bSJed Brown for (i=0; i<m; i++) { 249c4762a1bSJed Brown for (j=0; j<n; j++) { 250c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 251*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown } 254c4762a1bSJed Brown } 255c4762a1bSJed Brown } 256*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 257*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* a * dF/d(xdot) part */ 260c4762a1bSJed Brown if (adctx->Seed) 261c4762a1bSJed Brown fov_forward(tag2,m,n,p,u_vec,adctx->Seed,NULL,J); 262c4762a1bSJed Brown else 263c4762a1bSJed Brown jacobian(tag2,m,n,u_vec,J); 264c4762a1bSJed Brown if (adctx->sparse) { 265*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A,ADD_VALUES,m,p,adctx->Rec,J,&a)); 266c4762a1bSJed Brown } else { 267c4762a1bSJed Brown for (i=0; i<m; i++) { 268c4762a1bSJed Brown for (j=0; j<n; j++) { 269c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 270c4762a1bSJed Brown J[i][j] *= a; 271*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,1,&j,&J[i][j],ADD_VALUES)); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275c4762a1bSJed Brown } 276*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 277*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 278*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 279c4762a1bSJed Brown PetscFunctionReturn(0); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* 283c4762a1bSJed Brown Compute local portion of Jacobian for implicit TS in the special case where it is 284c4762a1bSJed Brown known that the mass matrix is simply the identity. i.e. We have 285c4762a1bSJed Brown a problem of the form 286c4762a1bSJed Brown du/dt = F(u). 287c4762a1bSJed Brown 288c4762a1bSJed Brown Input parameters: 289c4762a1bSJed Brown tag - tape identifier for df/dx part 290c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 291c4762a1bSJed Brown ctx - ADOL-C context, as defined above 292c4762a1bSJed Brown 293c4762a1bSJed Brown Output parameter: 294c4762a1bSJed Brown A - Mat object corresponding to Jacobian 295c4762a1bSJed Brown */ 296410585f6SHong Zhang PetscErrorCode PetscAdolcComputeIJacobianLocalIDMass(PetscInt tag,Mat A,const PetscScalar *u_vec,PetscReal a,void *ctx) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 299c4762a1bSJed Brown PetscInt i,j,m = adctx->m,n = adctx->n,p = adctx->p; 300c4762a1bSJed Brown PetscScalar **J; 301c4762a1bSJed Brown 302c4762a1bSJed Brown PetscFunctionBegin; 303*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* dF/dx part */ 306c4762a1bSJed Brown if (adctx->Seed) 307c4762a1bSJed Brown fov_forward(tag,m,n,p,u_vec,adctx->Seed,NULL,J); 308c4762a1bSJed Brown else 309c4762a1bSJed Brown jacobian(tag,m,n,u_vec,J); 310c4762a1bSJed Brown if (adctx->sparse) { 311*9566063dSJacob Faibussowitsch PetscCall(RecoverJacobianLocal(A,INSERT_VALUES,m,p,adctx->Rec,J,NULL)); 312c4762a1bSJed Brown } else { 313c4762a1bSJed Brown for (i=0; i<m; i++) { 314c4762a1bSJed Brown for (j=0; j<n; j++) { 315c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 316*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 317c4762a1bSJed Brown } 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } 320c4762a1bSJed Brown } 321*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 322*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 323*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* a * dF/d(xdot) part */ 326*9566063dSJacob Faibussowitsch PetscCall(MatShift(A,a)); 327c4762a1bSJed Brown PetscFunctionReturn(0); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* -------------------------------------------------------------------------------- 331c4762a1bSJed Brown Drivers for Jacobian w.r.t. a parameter 332c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* 335c4762a1bSJed Brown Compute Jacobian w.r.t a parameter for explicit TS. 336c4762a1bSJed Brown 337c4762a1bSJed Brown Input parameters: 338c4762a1bSJed Brown tag - tape identifier 339c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 340c4762a1bSJed Brown params - the parameters w.r.t. which we differentiate 341c4762a1bSJed Brown ctx - ADOL-C context, as defined above 342c4762a1bSJed Brown 343c4762a1bSJed Brown Output parameter: 344c4762a1bSJed Brown A - Mat object corresponding to Jacobian 345c4762a1bSJed Brown */ 346a8c08197SHong Zhang PetscErrorCode PetscAdolcComputeRHSJacobianP(PetscInt tag,Mat A,const PetscScalar *u_vec,PetscScalar *params,void *ctx) 347c4762a1bSJed Brown { 348c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 349c4762a1bSJed Brown PetscInt i,j = 0,m = adctx->m,n = adctx->n,p = adctx->num_params; 350c4762a1bSJed Brown PetscScalar **J,*concat,**S; 351c4762a1bSJed Brown 352c4762a1bSJed Brown PetscFunctionBegin; 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* Allocate memory and concatenate independent variable values with parameter */ 355*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 356*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n+p,&concat)); 357*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(n+p,p,&S)); 358*9566063dSJacob Faibussowitsch PetscCall(Subidentity(p,n,S)); 359c4762a1bSJed Brown for (i=0; i<n; i++) concat[i] = u_vec[i]; 360c4762a1bSJed Brown for (i=0; i<p; i++) concat[n+i] = params[i]; 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* Propagate the appropriate seed matrix through the forward mode of AD */ 363c4762a1bSJed Brown fov_forward(tag,m,n+p,p,concat,S,NULL,J); 364*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(S)); 365*9566063dSJacob Faibussowitsch PetscCall(PetscFree(concat)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* Set matrix values */ 368c4762a1bSJed Brown for (i=0; i<m; i++) { 369c4762a1bSJed Brown for (j=0; j<p; j++) { 370c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 371*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown } 374c4762a1bSJed Brown } 375*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 376*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 377*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 378c4762a1bSJed Brown PetscFunctionReturn(0); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 381c4762a1bSJed Brown /* 382c4762a1bSJed Brown Compute local portion of Jacobian w.r.t a parameter for explicit TS. 383c4762a1bSJed Brown 384c4762a1bSJed Brown Input parameters: 385c4762a1bSJed Brown tag - tape identifier 386c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 387c4762a1bSJed Brown params - the parameters w.r.t. which we differentiate 388c4762a1bSJed Brown ctx - ADOL-C context, as defined above 389c4762a1bSJed Brown 390c4762a1bSJed Brown Output parameter: 391c4762a1bSJed Brown A - Mat object corresponding to Jacobian 392c4762a1bSJed Brown */ 393a8c08197SHong Zhang PetscErrorCode PetscAdolcComputeRHSJacobianPLocal(PetscInt tag,Mat A,const PetscScalar *u_vec,PetscScalar *params,void *ctx) 394c4762a1bSJed Brown { 395c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 396c4762a1bSJed Brown PetscInt i,j = 0,m = adctx->m,n = adctx->n,p = adctx->num_params; 397c4762a1bSJed Brown PetscScalar **J,*concat,**S; 398c4762a1bSJed Brown 399c4762a1bSJed Brown PetscFunctionBegin; 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* Allocate memory and concatenate independent variable values with parameter */ 402*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 403*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n+p,&concat)); 404*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(n+p,p,&S)); 405*9566063dSJacob Faibussowitsch PetscCall(Subidentity(p,n,S)); 406c4762a1bSJed Brown for (i=0; i<n; i++) concat[i] = u_vec[i]; 407c4762a1bSJed Brown for (i=0; i<p; i++) concat[n+i] = params[i]; 408c4762a1bSJed Brown 409c4762a1bSJed Brown /* Propagate the appropriate seed matrix through the forward mode of AD */ 410c4762a1bSJed Brown fov_forward(tag,m,n+p,p,concat,S,NULL,J); 411*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(S)); 412*9566063dSJacob Faibussowitsch PetscCall(PetscFree(concat)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Set matrix values */ 415c4762a1bSJed Brown for (i=0; i<m; i++) { 416c4762a1bSJed Brown for (j=0; j<p; j++) { 417c4762a1bSJed Brown if (fabs(J[i][j]) > 1.e-16) { 418*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown } 422*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 423*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 424*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 425c4762a1bSJed Brown PetscFunctionReturn(0); 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* -------------------------------------------------------------------------------- 429c4762a1bSJed Brown Drivers for Jacobian diagonal 430c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* 433c4762a1bSJed Brown Compute local portion of Jacobian diagonal for implicit TS in compressed format and recover 434c4762a1bSJed Brown from this, using precomputed seed matrix and recovery vector. 435c4762a1bSJed Brown 436c4762a1bSJed Brown Input parameters: 437c4762a1bSJed Brown tag1 - tape identifier for df/dx part 438c4762a1bSJed Brown tag2 - tape identifier for df/d(xdot) part 439c4762a1bSJed Brown u_vec - vector at which to evaluate Jacobian 440c4762a1bSJed Brown ctx - ADOL-C context, as defined above 441c4762a1bSJed Brown 442c4762a1bSJed Brown Output parameter: 443c4762a1bSJed Brown diag - Vec object corresponding to Jacobian diagonal 444c4762a1bSJed Brown */ 445c4762a1bSJed Brown PetscErrorCode PetscAdolcComputeIJacobianAndDiagonalLocal(PetscInt tag1,PetscInt tag2,Vec diag,PetscScalar *u_vec,PetscReal a,void *ctx) 446c4762a1bSJed Brown { 447c4762a1bSJed Brown AdolcCtx *adctx = (AdolcCtx*)ctx; 448c4762a1bSJed Brown PetscInt i,m = adctx->m,n = adctx->n,p = adctx->p; 449c4762a1bSJed Brown PetscScalar **J; 450c4762a1bSJed Brown 451c4762a1bSJed Brown PetscFunctionBegin; 452*9566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(m,p,&J)); 453c4762a1bSJed Brown 454c4762a1bSJed Brown /* dF/dx part */ 455c4762a1bSJed Brown if (adctx->Seed) 456c4762a1bSJed Brown fov_forward(tag1,m,n,p,u_vec,adctx->Seed,NULL,J); 457c4762a1bSJed Brown else 458c4762a1bSJed Brown jacobian(tag1,m,n,u_vec,J); 459c4762a1bSJed Brown if (adctx->sparse) { 460*9566063dSJacob Faibussowitsch PetscCall(RecoverDiagonalLocal(diag,INSERT_VALUES,m,adctx->rec,J,NULL)); 461c4762a1bSJed Brown } else { 462c4762a1bSJed Brown for (i=0; i<m; i++) { 463c4762a1bSJed Brown if (fabs(J[i][i]) > 1.e-16) { 464*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(diag,1,&i,&J[i][i],INSERT_VALUES)); 465c4762a1bSJed Brown } 466c4762a1bSJed Brown } 467c4762a1bSJed Brown } 468*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(diag)); 469*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(diag)); 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* a * dF/d(xdot) part */ 472c4762a1bSJed Brown if (adctx->Seed) 473c4762a1bSJed Brown fov_forward(tag2,m,n,p,u_vec,adctx->Seed,NULL,J); 474c4762a1bSJed Brown else 475c4762a1bSJed Brown jacobian(tag2,m,n,u_vec,J); 476c4762a1bSJed Brown if (adctx->sparse) { 477*9566063dSJacob Faibussowitsch PetscCall(RecoverDiagonalLocal(diag,ADD_VALUES,m,adctx->rec,J,NULL)); 478c4762a1bSJed Brown } else { 479c4762a1bSJed Brown for (i=0; i<m; i++) { 480c4762a1bSJed Brown if (fabs(J[i][i]) > 1.e-16) { 481c4762a1bSJed Brown J[i][i] *= a; 482*9566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(diag,1,&i,&J[i][i],ADD_VALUES)); 483c4762a1bSJed Brown } 484c4762a1bSJed Brown } 485c4762a1bSJed Brown } 486*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(diag)); 487*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(diag)); 488*9566063dSJacob Faibussowitsch PetscCall(AdolcFree2(J)); 489c4762a1bSJed Brown PetscFunctionReturn(0); 490c4762a1bSJed Brown } 491