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