1c4762a1bSJed Brown #include <petscdmda.h>
2c4762a1bSJed Brown #include <petscts.h>
3c4762a1bSJed Brown #include <adolc/interfaces.h>
4c4762a1bSJed Brown #include "contexts.cxx"
5c4762a1bSJed Brown
6c4762a1bSJed Brown /*
7c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc.
8c4762a1bSJed Brown
9c4762a1bSJed Brown For documentation on ADOL-C, see
10c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
11c4762a1bSJed Brown */
12c4762a1bSJed Brown
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown ADOL-C implementation for Jacobian vector product, using the forward mode of AD.
15c4762a1bSJed Brown Intended to overload MatMult in matrix-free methods where implicit timestepping
16c4762a1bSJed Brown has been used.
17c4762a1bSJed Brown
18c4762a1bSJed Brown For an implicit Jacobian we may use the rule that
19c4762a1bSJed Brown G = M*xdot + f(x) ==> dG/dx = a*M + df/dx,
20c4762a1bSJed Brown where a = d(xdot)/dx is a constant. Evaluated at x0 and acting upon a vector x1:
21c4762a1bSJed Brown (dG/dx)(x0) * x1 = (a*df/d(xdot)(x0) + df/dx)(x0) * x1.
22c4762a1bSJed Brown
23c4762a1bSJed Brown Input parameters:
24c4762a1bSJed Brown A_shell - Jacobian matrix of MatShell type
25c4762a1bSJed Brown X - vector to be multiplied by A_shell
26c4762a1bSJed Brown
27c4762a1bSJed Brown Output parameters:
28c4762a1bSJed Brown Y - product of A_shell and X
29c4762a1bSJed Brown */
PetscAdolcIJacobianVectorProduct(Mat A_shell,Vec X,Vec Y)30d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianVectorProduct(Mat A_shell, Vec X, Vec Y)
31d71ae5a4SJacob Faibussowitsch {
32c4762a1bSJed Brown AdolcMatCtx *mctx;
33c4762a1bSJed Brown PetscInt m, n, i, j, k = 0, d;
34c4762a1bSJed Brown const PetscScalar *x0;
35c4762a1bSJed Brown PetscScalar *action, *x1;
36c4762a1bSJed Brown Vec localX1;
37c4762a1bSJed Brown DM da;
38c4762a1bSJed Brown DMDALocalInfo info;
39c4762a1bSJed Brown
40c4762a1bSJed Brown PetscFunctionBegin;
41c4762a1bSJed Brown /* Get matrix-free context info */
429566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A_shell, &mctx));
43c4762a1bSJed Brown m = mctx->m;
44c4762a1bSJed Brown n = mctx->n;
45c4762a1bSJed Brown
46c4762a1bSJed Brown /* Get local input vectors and extract data, x0 and x1*/
479566063dSJacob Faibussowitsch PetscCall(TSGetDM(mctx->ts, &da));
489566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
499566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX1));
509566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX1));
519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX1));
52c4762a1bSJed Brown
539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mctx->localX0, &x0));
549566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX1, &x1));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* dF/dx part */
579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &action));
589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event1, 0, 0, 0, 0));
59c4762a1bSJed Brown fos_forward(mctx->tag1, m, n, 0, x0, x1, NULL, action);
60c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
61c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
62c4762a1bSJed Brown for (d = 0; d < 2; d++) {
6348a46eb9SPierre Jolivet if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) PetscCall(VecSetValuesLocal(Y, 1, &k, &action[k], INSERT_VALUES));
64c4762a1bSJed Brown k++;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown }
67c4762a1bSJed Brown }
689566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event1, 0, 0, 0, 0));
69c4762a1bSJed Brown k = 0;
709566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
719566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y)); /* to INSERT_VALUES and ADD_VALUES */
72c4762a1bSJed Brown
73c4762a1bSJed Brown /* a * dF/d(xdot) part */
749566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event2, 0, 0, 0, 0));
75c4762a1bSJed Brown fos_forward(mctx->tag2, m, n, 0, x0, x1, NULL, action);
76c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
77c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
78c4762a1bSJed Brown for (d = 0; d < 2; d++) {
79c4762a1bSJed Brown if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) {
80c4762a1bSJed Brown action[k] *= mctx->shift;
819566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(Y, 1, &k, &action[k], ADD_VALUES));
82c4762a1bSJed Brown }
83c4762a1bSJed Brown k++;
84c4762a1bSJed Brown }
85c4762a1bSJed Brown }
86c4762a1bSJed Brown }
879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event2, 0, 0, 0, 0));
889566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y));
899566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y));
909566063dSJacob Faibussowitsch PetscCall(PetscFree(action));
91c4762a1bSJed Brown
92c4762a1bSJed Brown /* Restore local vector */
939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX1, &x1));
949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mctx->localX0, &x0));
959566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX1));
96*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown ADOL-C implementation for Jacobian vector product, using the forward mode of AD.
101c4762a1bSJed Brown Intended to overload MatMult in matrix-free methods where implicit timestepping
102c4762a1bSJed Brown has been applied to a problem of the form
103c4762a1bSJed Brown du/dt = F(u).
104c4762a1bSJed Brown
105c4762a1bSJed Brown Input parameters:
106c4762a1bSJed Brown A_shell - Jacobian matrix of MatShell type
107c4762a1bSJed Brown X - vector to be multiplied by A_shell
108c4762a1bSJed Brown
109c4762a1bSJed Brown Output parameters:
110c4762a1bSJed Brown Y - product of A_shell and X
111c4762a1bSJed Brown */
PetscAdolcIJacobianVectorProductIDMass(Mat A_shell,Vec X,Vec Y)112d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianVectorProductIDMass(Mat A_shell, Vec X, Vec Y)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown AdolcMatCtx *mctx;
115c4762a1bSJed Brown PetscInt m, n, i, j, k = 0, d;
116c4762a1bSJed Brown const PetscScalar *x0;
117c4762a1bSJed Brown PetscScalar *action, *x1;
118c4762a1bSJed Brown Vec localX1;
119c4762a1bSJed Brown DM da;
120c4762a1bSJed Brown DMDALocalInfo info;
121c4762a1bSJed Brown
122c4762a1bSJed Brown PetscFunctionBegin;
123c4762a1bSJed Brown /* Get matrix-free context info */
1249566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A_shell, &mctx));
125c4762a1bSJed Brown m = mctx->m;
126c4762a1bSJed Brown n = mctx->n;
127c4762a1bSJed Brown
128c4762a1bSJed Brown /* Get local input vectors and extract data, x0 and x1*/
1299566063dSJacob Faibussowitsch PetscCall(TSGetDM(mctx->ts, &da));
1309566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
1319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX1));
1329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX1));
1339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX1));
134c4762a1bSJed Brown
1359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mctx->localX0, &x0));
1369566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX1, &x1));
137c4762a1bSJed Brown
138c4762a1bSJed Brown /* dF/dx part */
1399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &action));
1409566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event1, 0, 0, 0, 0));
141c4762a1bSJed Brown fos_forward(mctx->tag1, m, n, 0, x0, x1, NULL, action);
142c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
143c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
144c4762a1bSJed Brown for (d = 0; d < 2; d++) {
14548a46eb9SPierre Jolivet if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) PetscCall(VecSetValuesLocal(Y, 1, &k, &action[k], INSERT_VALUES));
146c4762a1bSJed Brown k++;
147c4762a1bSJed Brown }
148c4762a1bSJed Brown }
149c4762a1bSJed Brown }
1509566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event1, 0, 0, 0, 0));
151c4762a1bSJed Brown k = 0;
1529566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
1539566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y)); /* to INSERT_VALUES and ADD_VALUES */
1549566063dSJacob Faibussowitsch PetscCall(PetscFree(action));
155c4762a1bSJed Brown
156c4762a1bSJed Brown /* Restore local vector */
1579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX1, &x1));
1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mctx->localX0, &x0));
1599566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX1));
160c4762a1bSJed Brown
161c4762a1bSJed Brown /* a * dF/d(xdot) part */
1629566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event2, 0, 0, 0, 0));
1639566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, mctx->shift, X));
1649566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event2, 0, 0, 0, 0));
165*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown
168c4762a1bSJed Brown /*
169c4762a1bSJed Brown ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
170c4762a1bSJed Brown Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
171c4762a1bSJed Brown has been used.
172c4762a1bSJed Brown
173c4762a1bSJed Brown Input parameters:
174c4762a1bSJed Brown A_shell - Jacobian matrix of MatShell type
175c4762a1bSJed Brown Y - vector to be multiplied by A_shell transpose
176c4762a1bSJed Brown
177c4762a1bSJed Brown Output parameters:
178c4762a1bSJed Brown X - product of A_shell transpose and X
179c4762a1bSJed Brown */
PetscAdolcIJacobianTransposeVectorProduct(Mat A_shell,Vec Y,Vec X)180d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianTransposeVectorProduct(Mat A_shell, Vec Y, Vec X)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown AdolcMatCtx *mctx;
183c4762a1bSJed Brown PetscInt m, n, i, j, k = 0, d;
184c4762a1bSJed Brown const PetscScalar *x;
185c4762a1bSJed Brown PetscScalar *action, *y;
186c4762a1bSJed Brown Vec localY;
187c4762a1bSJed Brown DM da;
188c4762a1bSJed Brown DMDALocalInfo info;
189c4762a1bSJed Brown
190c4762a1bSJed Brown PetscFunctionBegin;
191c4762a1bSJed Brown /* Get matrix-free context info */
1929566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A_shell, &mctx));
193c4762a1bSJed Brown m = mctx->m;
194c4762a1bSJed Brown n = mctx->n;
195c4762a1bSJed Brown
196c4762a1bSJed Brown /* Get local input vectors and extract data, x0 and x1*/
1979566063dSJacob Faibussowitsch PetscCall(TSGetDM(mctx->ts, &da));
1989566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
1999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localY));
2009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, Y, INSERT_VALUES, localY));
2019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, Y, INSERT_VALUES, localY));
2029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mctx->localX0, &x));
2039566063dSJacob Faibussowitsch PetscCall(VecGetArray(localY, &y));
204c4762a1bSJed Brown
205c4762a1bSJed Brown /* dF/dx part */
2069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &action));
2079566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event3, 0, 0, 0, 0));
2089371c9d4SSatish Balay if (!mctx->flg) zos_forward(mctx->tag1, m, n, 1, x, NULL);
209c4762a1bSJed Brown fos_reverse(mctx->tag1, m, n, y, action);
210c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
211c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
212c4762a1bSJed Brown for (d = 0; d < 2; d++) {
21348a46eb9SPierre Jolivet if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) PetscCall(VecSetValuesLocal(X, 1, &k, &action[k], INSERT_VALUES));
214c4762a1bSJed Brown k++;
215c4762a1bSJed Brown }
216c4762a1bSJed Brown }
217c4762a1bSJed Brown }
2189566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event3, 0, 0, 0, 0));
219c4762a1bSJed Brown k = 0;
2209566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
2219566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); /* to INSERT_VALUES and ADD_VALUES */
222c4762a1bSJed Brown
223c4762a1bSJed Brown /* a * dF/d(xdot) part */
2249566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event4, 0, 0, 0, 0));
225c4762a1bSJed Brown if (!mctx->flg) {
226c4762a1bSJed Brown zos_forward(mctx->tag2, m, n, 1, x, NULL);
227c4762a1bSJed Brown mctx->flg = PETSC_TRUE;
228c4762a1bSJed Brown }
229c4762a1bSJed Brown fos_reverse(mctx->tag2, m, n, y, action);
230c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
231c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
232c4762a1bSJed Brown for (d = 0; d < 2; d++) {
233c4762a1bSJed Brown if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) {
234c4762a1bSJed Brown action[k] *= mctx->shift;
2359566063dSJacob Faibussowitsch PetscCall(VecSetValuesLocal(X, 1, &k, &action[k], ADD_VALUES));
236c4762a1bSJed Brown }
237c4762a1bSJed Brown k++;
238c4762a1bSJed Brown }
239c4762a1bSJed Brown }
240c4762a1bSJed Brown }
2419566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event4, 0, 0, 0, 0));
2429566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X));
2439566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X));
2449566063dSJacob Faibussowitsch PetscCall(PetscFree(action));
245c4762a1bSJed Brown
246c4762a1bSJed Brown /* Restore local vector */
2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localY, &y));
2489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mctx->localX0, &x));
2499566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localY));
250*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
251c4762a1bSJed Brown }
252c4762a1bSJed Brown
253c4762a1bSJed Brown /*
254c4762a1bSJed Brown ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
255c4762a1bSJed Brown Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
256c4762a1bSJed Brown has been applied to a problem of the form
257c4762a1bSJed Brown du/dt = F(u).
258c4762a1bSJed Brown
259c4762a1bSJed Brown Input parameters:
260c4762a1bSJed Brown A_shell - Jacobian matrix of MatShell type
261c4762a1bSJed Brown Y - vector to be multiplied by A_shell transpose
262c4762a1bSJed Brown
263c4762a1bSJed Brown Output parameters:
264c4762a1bSJed Brown X - product of A_shell transpose and X
265c4762a1bSJed Brown */
PetscAdolcIJacobianTransposeVectorProductIDMass(Mat A_shell,Vec Y,Vec X)266d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianTransposeVectorProductIDMass(Mat A_shell, Vec Y, Vec X)
267d71ae5a4SJacob Faibussowitsch {
268c4762a1bSJed Brown AdolcMatCtx *mctx;
269c4762a1bSJed Brown PetscInt m, n, i, j, k = 0, d;
270c4762a1bSJed Brown const PetscScalar *x;
271c4762a1bSJed Brown PetscScalar *action, *y;
272c4762a1bSJed Brown Vec localY;
273c4762a1bSJed Brown DM da;
274c4762a1bSJed Brown DMDALocalInfo info;
275c4762a1bSJed Brown
276c4762a1bSJed Brown PetscFunctionBegin;
277c4762a1bSJed Brown /* Get matrix-free context info */
2789566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A_shell, &mctx));
279c4762a1bSJed Brown m = mctx->m;
280c4762a1bSJed Brown n = mctx->n;
281c4762a1bSJed Brown
282c4762a1bSJed Brown /* Get local input vectors and extract data, x0 and x1*/
2839566063dSJacob Faibussowitsch PetscCall(TSGetDM(mctx->ts, &da));
2849566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
2859566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localY));
2869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, Y, INSERT_VALUES, localY));
2879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, Y, INSERT_VALUES, localY));
2889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mctx->localX0, &x));
2899566063dSJacob Faibussowitsch PetscCall(VecGetArray(localY, &y));
290c4762a1bSJed Brown
291c4762a1bSJed Brown /* dF/dx part */
2929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &action));
2939566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event3, 0, 0, 0, 0));
294c4762a1bSJed Brown if (!mctx->flg) zos_forward(mctx->tag1, m, n, 1, x, NULL);
295c4762a1bSJed Brown fos_reverse(mctx->tag1, m, n, y, action);
296c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) {
297c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) {
298c4762a1bSJed Brown for (d = 0; d < 2; d++) {
29948a46eb9SPierre Jolivet if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) PetscCall(VecSetValuesLocal(X, 1, &k, &action[k], INSERT_VALUES));
300c4762a1bSJed Brown k++;
301c4762a1bSJed Brown }
302c4762a1bSJed Brown }
303c4762a1bSJed Brown }
3049566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event3, 0, 0, 0, 0));
305c4762a1bSJed Brown k = 0;
3069566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
3079566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); /* to INSERT_VALUES and ADD_VALUES */
3089566063dSJacob Faibussowitsch PetscCall(PetscFree(action));
309c4762a1bSJed Brown
310c4762a1bSJed Brown /* Restore local vector */
3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localY, &y));
3129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mctx->localX0, &x));
3139566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localY));
314c4762a1bSJed Brown
315c4762a1bSJed Brown /* a * dF/d(xdot) part */
3169566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(mctx->event4, 0, 0, 0, 0));
3179566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, mctx->shift, Y));
3189566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(mctx->event4, 0, 0, 0, 0));
319*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
320c4762a1bSJed Brown }
321