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