xref: /petsc/src/ts/tutorials/autodiff/adolc-utils/drivers.cxx (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 */
31*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeRHSJacobian(PetscInt tag, Mat A, const PetscScalar *u_vec, void *ctx)
32*d71ae5a4SJacob Faibussowitsch {
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;
389566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
399371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J);
409371c9d4SSatish Balay   else jacobian(tag, m, n, u_vec, J);
41c4762a1bSJed Brown   if (adctx->sparse) {
429566063dSJacob Faibussowitsch     PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
43c4762a1bSJed Brown   } else {
44c4762a1bSJed Brown     for (i = 0; i < m; i++) {
45c4762a1bSJed Brown       for (j = 0; j < n; j++) {
4648a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
47c4762a1bSJed Brown       }
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown   }
509566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown   Compute Jacobian for explicit TS in compressed format and recover from this, using
58c4762a1bSJed Brown   precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is
59c4762a1bSJed Brown   assembled (not recommended for non-toy problems!).
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   Input parameters:
62c4762a1bSJed Brown   tag   - tape identifier
63c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
64c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   Output parameter:
67c4762a1bSJed Brown   A     - Mat object corresponding to Jacobian
68c4762a1bSJed Brown */
69*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeRHSJacobianLocal(PetscInt tag, Mat A, const PetscScalar *u_vec, void *ctx)
70*d71ae5a4SJacob Faibussowitsch {
71c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
72c4762a1bSJed Brown   PetscInt      i, j, m = adctx->m, n = adctx->n, p = adctx->p;
73c4762a1bSJed Brown   PetscScalar **J;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
779371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J);
789371c9d4SSatish Balay   else jacobian(tag, m, n, u_vec, J);
79c4762a1bSJed Brown   if (adctx->sparse) {
809566063dSJacob Faibussowitsch     PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
81c4762a1bSJed Brown   } else {
82c4762a1bSJed Brown     for (i = 0; i < m; i++) {
83c4762a1bSJed Brown       for (j = 0; j < n; j++) {
8448a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
85c4762a1bSJed Brown       }
86c4762a1bSJed Brown     }
87c4762a1bSJed Brown   }
889566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
91c4762a1bSJed Brown   PetscFunctionReturn(0);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown /*
95c4762a1bSJed Brown   Compute Jacobian for implicit TS in compressed format and recover from this, using
96c4762a1bSJed Brown   precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is
97c4762a1bSJed Brown   assembled (not recommended for non-toy problems!).
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   Input parameters:
100c4762a1bSJed Brown   tag1   - tape identifier for df/dx part
101c4762a1bSJed Brown   tag2   - tape identifier for df/d(xdot) part
102c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
103c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   Output parameter:
106c4762a1bSJed Brown   A     - Mat object corresponding to Jacobian
107c4762a1bSJed Brown */
108*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeIJacobian(PetscInt tag1, PetscInt tag2, Mat A, const PetscScalar *u_vec, PetscReal a, void *ctx)
109*d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
111c4762a1bSJed Brown   PetscInt      i, j, m = adctx->m, n = adctx->n, p = adctx->p;
112c4762a1bSJed Brown   PetscScalar **J;
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* dF/dx part */
1189371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J);
1199371c9d4SSatish Balay   else jacobian(tag1, m, n, u_vec, J);
1209566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(A));
121c4762a1bSJed Brown   if (adctx->sparse) {
1229566063dSJacob Faibussowitsch     PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
123c4762a1bSJed Brown   } else {
124c4762a1bSJed Brown     for (i = 0; i < m; i++) {
125c4762a1bSJed Brown       for (j = 0; j < n; j++) {
12648a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
127c4762a1bSJed Brown       }
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown   }
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* a * dF/d(xdot) part */
1349371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J);
1359371c9d4SSatish Balay   else jacobian(tag2, m, n, u_vec, J);
136c4762a1bSJed Brown   if (adctx->sparse) {
1379566063dSJacob Faibussowitsch     PetscCall(RecoverJacobian(A, ADD_VALUES, m, p, adctx->Rec, J, &a));
138c4762a1bSJed Brown   } else {
139c4762a1bSJed Brown     for (i = 0; i < m; i++) {
140c4762a1bSJed Brown       for (j = 0; j < n; j++) {
141c4762a1bSJed Brown         if (fabs(J[i][j]) > 1.e-16) {
142c4762a1bSJed Brown           J[i][j] *= a;
1439566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], ADD_VALUES));
144c4762a1bSJed Brown         }
145c4762a1bSJed Brown       }
146c4762a1bSJed Brown     }
147c4762a1bSJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
151c4762a1bSJed Brown   PetscFunctionReturn(0);
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
154c4762a1bSJed Brown /*
155c4762a1bSJed Brown   Compute Jacobian for implicit TS in the special case where it is
156c4762a1bSJed Brown   known that the mass matrix is simply the identity. i.e. We have
157c4762a1bSJed Brown   a problem of the form
158c4762a1bSJed Brown                         du/dt = F(u).
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   Input parameters:
161c4762a1bSJed Brown   tag   - tape identifier for df/dx part
162c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
163c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   Output parameter:
166c4762a1bSJed Brown   A     - Mat object corresponding to Jacobian
167c4762a1bSJed Brown */
168*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeIJacobianIDMass(PetscInt tag, Mat A, PetscScalar *u_vec, PetscReal a, void *ctx)
169*d71ae5a4SJacob Faibussowitsch {
170c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
171c4762a1bSJed Brown   PetscInt      i, j, m = adctx->m, n = adctx->n, p = adctx->p;
172c4762a1bSJed Brown   PetscScalar **J;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* dF/dx part */
1789371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J);
1799371c9d4SSatish Balay   else jacobian(tag, m, n, u_vec, J);
1809566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(A));
181c4762a1bSJed Brown   if (adctx->sparse) {
1829566063dSJacob Faibussowitsch     PetscCall(RecoverJacobian(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
183c4762a1bSJed Brown   } else {
184c4762a1bSJed Brown     for (i = 0; i < m; i++) {
185c4762a1bSJed Brown       for (j = 0; j < n; j++) {
18648a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
187c4762a1bSJed Brown       }
188c4762a1bSJed Brown     }
189c4762a1bSJed Brown   }
1909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1929566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* a * dF/d(xdot) part */
1959566063dSJacob Faibussowitsch   PetscCall(MatShift(A, a));
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /*
200c4762a1bSJed Brown   Compute local portion of Jacobian for implicit TS in compressed format and recover from this, using
201c4762a1bSJed Brown   precomputed seed and recovery matrices. If sparse mode is not used, full Jacobian is
202c4762a1bSJed Brown   assembled (not recommended for non-toy problems!).
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   Input parameters:
205c4762a1bSJed Brown   tag1   - tape identifier for df/dx part
206c4762a1bSJed Brown   tag2   - tape identifier for df/d(xdot) part
207c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
208c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   Output parameter:
211c4762a1bSJed Brown   A     - Mat object corresponding to Jacobian
212c4762a1bSJed Brown */
213*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeIJacobianLocal(PetscInt tag1, PetscInt tag2, Mat A, PetscScalar *u_vec, PetscReal a, void *ctx)
214*d71ae5a4SJacob Faibussowitsch {
215c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
216c4762a1bSJed Brown   PetscInt      i, j, m = adctx->m, n = adctx->n, p = adctx->p;
217c4762a1bSJed Brown   PetscScalar **J;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /* dF/dx part */
2239371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J);
2249371c9d4SSatish Balay   else jacobian(tag1, m, n, u_vec, J);
225c4762a1bSJed Brown   if (adctx->sparse) {
2269566063dSJacob Faibussowitsch     PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
227c4762a1bSJed Brown   } else {
228c4762a1bSJed Brown     for (i = 0; i < m; i++) {
229c4762a1bSJed Brown       for (j = 0; j < n; j++) {
23048a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
231c4762a1bSJed Brown       }
232c4762a1bSJed Brown     }
233c4762a1bSJed Brown   }
2349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* a * dF/d(xdot) part */
2389371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J);
2399371c9d4SSatish Balay   else jacobian(tag2, m, n, u_vec, J);
240c4762a1bSJed Brown   if (adctx->sparse) {
2419566063dSJacob Faibussowitsch     PetscCall(RecoverJacobianLocal(A, ADD_VALUES, m, p, adctx->Rec, J, &a));
242c4762a1bSJed Brown   } else {
243c4762a1bSJed Brown     for (i = 0; i < m; i++) {
244c4762a1bSJed Brown       for (j = 0; j < n; j++) {
245c4762a1bSJed Brown         if (fabs(J[i][j]) > 1.e-16) {
246c4762a1bSJed Brown           J[i][j] *= a;
2479566063dSJacob Faibussowitsch           PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], ADD_VALUES));
248c4762a1bSJed Brown         }
249c4762a1bSJed Brown       }
250c4762a1bSJed Brown     }
251c4762a1bSJed Brown   }
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
255c4762a1bSJed Brown   PetscFunctionReturn(0);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /*
259c4762a1bSJed Brown   Compute local portion of Jacobian for implicit TS in the special case where it is
260c4762a1bSJed Brown   known that the mass matrix is simply the identity. i.e. We have
261c4762a1bSJed Brown   a problem of the form
262c4762a1bSJed Brown                         du/dt = F(u).
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   Input parameters:
265c4762a1bSJed Brown   tag   - tape identifier for df/dx part
266c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
267c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   Output parameter:
270c4762a1bSJed Brown   A     - Mat object corresponding to Jacobian
271c4762a1bSJed Brown */
272*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeIJacobianLocalIDMass(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscReal a, void *ctx)
273*d71ae5a4SJacob Faibussowitsch {
274c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
275c4762a1bSJed Brown   PetscInt      i, j, m = adctx->m, n = adctx->n, p = adctx->p;
276c4762a1bSJed Brown   PetscScalar **J;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   PetscFunctionBegin;
2799566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /* dF/dx part */
2829371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag, m, n, p, u_vec, adctx->Seed, NULL, J);
2839371c9d4SSatish Balay   else jacobian(tag, m, n, u_vec, J);
284c4762a1bSJed Brown   if (adctx->sparse) {
2859566063dSJacob Faibussowitsch     PetscCall(RecoverJacobianLocal(A, INSERT_VALUES, m, p, adctx->Rec, J, NULL));
286c4762a1bSJed Brown   } else {
287c4762a1bSJed Brown     for (i = 0; i < m; i++) {
288c4762a1bSJed Brown       for (j = 0; j < n; j++) {
28948a46eb9SPierre Jolivet         if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
290c4762a1bSJed Brown       }
291c4762a1bSJed Brown     }
292c4762a1bSJed Brown   }
2939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2959566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /* a * dF/d(xdot) part */
2989566063dSJacob Faibussowitsch   PetscCall(MatShift(A, a));
299c4762a1bSJed Brown   PetscFunctionReturn(0);
300c4762a1bSJed Brown }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown /* --------------------------------------------------------------------------------
303c4762a1bSJed Brown    Drivers for Jacobian w.r.t. a parameter
304c4762a1bSJed Brown    ----------------------------------------------------------------------------- */
305c4762a1bSJed Brown 
306c4762a1bSJed Brown /*
307c4762a1bSJed Brown   Compute Jacobian w.r.t a parameter for explicit TS.
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   Input parameters:
310c4762a1bSJed Brown   tag    - tape identifier
311c4762a1bSJed Brown   u_vec  - vector at which to evaluate Jacobian
312c4762a1bSJed Brown   params - the parameters w.r.t. which we differentiate
313c4762a1bSJed Brown   ctx    - ADOL-C context, as defined above
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   Output parameter:
316c4762a1bSJed Brown   A      - Mat object corresponding to Jacobian
317c4762a1bSJed Brown */
318*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeRHSJacobianP(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscScalar *params, void *ctx)
319*d71ae5a4SJacob Faibussowitsch {
320c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
321c4762a1bSJed Brown   PetscInt      i, j = 0, m = adctx->m, n = adctx->n, p = adctx->num_params;
322c4762a1bSJed Brown   PetscScalar **J, *concat, **S;
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   PetscFunctionBegin;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Allocate memory and concatenate independent variable values with parameter */
3279566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
3289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + p, &concat));
3299566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(n + p, p, &S));
3309566063dSJacob Faibussowitsch   PetscCall(Subidentity(p, n, S));
331c4762a1bSJed Brown   for (i = 0; i < n; i++) concat[i] = u_vec[i];
332c4762a1bSJed Brown   for (i = 0; i < p; i++) concat[n + i] = params[i];
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* Propagate the appropriate seed matrix through the forward mode of AD */
335c4762a1bSJed Brown   fov_forward(tag, m, n + p, p, concat, S, NULL, J);
3369566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(S));
3379566063dSJacob Faibussowitsch   PetscCall(PetscFree(concat));
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   /* Set matrix values */
340c4762a1bSJed Brown   for (i = 0; i < m; i++) {
341c4762a1bSJed Brown     for (j = 0; j < p; j++) {
34248a46eb9SPierre Jolivet       if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
343c4762a1bSJed Brown     }
344c4762a1bSJed Brown   }
3459566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
3469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
348c4762a1bSJed Brown   PetscFunctionReturn(0);
349c4762a1bSJed Brown }
350c4762a1bSJed Brown 
351c4762a1bSJed Brown /*
352c4762a1bSJed Brown   Compute local portion of Jacobian w.r.t a parameter for explicit TS.
353c4762a1bSJed Brown 
354c4762a1bSJed Brown   Input parameters:
355c4762a1bSJed Brown   tag    - tape identifier
356c4762a1bSJed Brown   u_vec  - vector at which to evaluate Jacobian
357c4762a1bSJed Brown   params - the parameters w.r.t. which we differentiate
358c4762a1bSJed Brown   ctx    - ADOL-C context, as defined above
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   Output parameter:
361c4762a1bSJed Brown   A      - Mat object corresponding to Jacobian
362c4762a1bSJed Brown */
363*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeRHSJacobianPLocal(PetscInt tag, Mat A, const PetscScalar *u_vec, PetscScalar *params, void *ctx)
364*d71ae5a4SJacob Faibussowitsch {
365c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
366c4762a1bSJed Brown   PetscInt      i, j = 0, m = adctx->m, n = adctx->n, p = adctx->num_params;
367c4762a1bSJed Brown   PetscScalar **J, *concat, **S;
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   PetscFunctionBegin;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /* Allocate memory and concatenate independent variable values with parameter */
3729566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
3739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + p, &concat));
3749566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(n + p, p, &S));
3759566063dSJacob Faibussowitsch   PetscCall(Subidentity(p, n, S));
376c4762a1bSJed Brown   for (i = 0; i < n; i++) concat[i] = u_vec[i];
377c4762a1bSJed Brown   for (i = 0; i < p; i++) concat[n + i] = params[i];
378c4762a1bSJed Brown 
379c4762a1bSJed Brown   /* Propagate the appropriate seed matrix through the forward mode of AD */
380c4762a1bSJed Brown   fov_forward(tag, m, n + p, p, concat, S, NULL, J);
3819566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(S));
3829566063dSJacob Faibussowitsch   PetscCall(PetscFree(concat));
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   /* Set matrix values */
385c4762a1bSJed Brown   for (i = 0; i < m; i++) {
386c4762a1bSJed Brown     for (j = 0; j < p; j++) {
38748a46eb9SPierre Jolivet       if (fabs(J[i][j]) > 1.e-16) PetscCall(MatSetValuesLocal(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES));
388c4762a1bSJed Brown     }
389c4762a1bSJed Brown   }
3909566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
3919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
393c4762a1bSJed Brown   PetscFunctionReturn(0);
394c4762a1bSJed Brown }
395c4762a1bSJed Brown 
396c4762a1bSJed Brown /* --------------------------------------------------------------------------------
397c4762a1bSJed Brown    Drivers for Jacobian diagonal
398c4762a1bSJed Brown    ----------------------------------------------------------------------------- */
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /*
401c4762a1bSJed Brown   Compute local portion of Jacobian diagonal for implicit TS in compressed format and recover
402c4762a1bSJed Brown   from this, using precomputed seed matrix and recovery vector.
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   Input parameters:
405c4762a1bSJed Brown   tag1  - tape identifier for df/dx part
406c4762a1bSJed Brown   tag2  - tape identifier for df/d(xdot) part
407c4762a1bSJed Brown   u_vec - vector at which to evaluate Jacobian
408c4762a1bSJed Brown   ctx   - ADOL-C context, as defined above
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   Output parameter:
411c4762a1bSJed Brown   diag  - Vec object corresponding to Jacobian diagonal
412c4762a1bSJed Brown */
413*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcComputeIJacobianAndDiagonalLocal(PetscInt tag1, PetscInt tag2, Vec diag, PetscScalar *u_vec, PetscReal a, void *ctx)
414*d71ae5a4SJacob Faibussowitsch {
415c4762a1bSJed Brown   AdolcCtx     *adctx = (AdolcCtx *)ctx;
416c4762a1bSJed Brown   PetscInt      i, m = adctx->m, n = adctx->n, p = adctx->p;
417c4762a1bSJed Brown   PetscScalar **J;
418c4762a1bSJed Brown 
419c4762a1bSJed Brown   PetscFunctionBegin;
4209566063dSJacob Faibussowitsch   PetscCall(AdolcMalloc2(m, p, &J));
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   /* dF/dx part */
4239371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag1, m, n, p, u_vec, adctx->Seed, NULL, J);
4249371c9d4SSatish Balay   else jacobian(tag1, m, n, u_vec, J);
425c4762a1bSJed Brown   if (adctx->sparse) {
4269566063dSJacob Faibussowitsch     PetscCall(RecoverDiagonalLocal(diag, INSERT_VALUES, m, adctx->rec, J, NULL));
427c4762a1bSJed Brown   } else {
428c4762a1bSJed Brown     for (i = 0; i < m; i++) {
42948a46eb9SPierre Jolivet       if (fabs(J[i][i]) > 1.e-16) PetscCall(VecSetValuesLocal(diag, 1, &i, &J[i][i], INSERT_VALUES));
430c4762a1bSJed Brown     }
431c4762a1bSJed Brown   }
4329566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(diag));
4339566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(diag));
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   /* a * dF/d(xdot) part */
4369371c9d4SSatish Balay   if (adctx->Seed) fov_forward(tag2, m, n, p, u_vec, adctx->Seed, NULL, J);
4379371c9d4SSatish Balay   else jacobian(tag2, m, n, u_vec, J);
438c4762a1bSJed Brown   if (adctx->sparse) {
4399566063dSJacob Faibussowitsch     PetscCall(RecoverDiagonalLocal(diag, ADD_VALUES, m, adctx->rec, J, NULL));
440c4762a1bSJed Brown   } else {
441c4762a1bSJed Brown     for (i = 0; i < m; i++) {
442c4762a1bSJed Brown       if (fabs(J[i][i]) > 1.e-16) {
443c4762a1bSJed Brown         J[i][i] *= a;
4449566063dSJacob Faibussowitsch         PetscCall(VecSetValuesLocal(diag, 1, &i, &J[i][i], ADD_VALUES));
445c4762a1bSJed Brown       }
446c4762a1bSJed Brown     }
447c4762a1bSJed Brown   }
4489566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(diag));
4499566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(diag));
4509566063dSJacob Faibussowitsch   PetscCall(AdolcFree2(J));
451c4762a1bSJed Brown   PetscFunctionReturn(0);
452c4762a1bSJed Brown }
453