xref: /petsc/src/ts/tutorials/autodiff/adolc-utils/matfree.cxx (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 */
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 
42c4762a1bSJed Brown   /* Get matrix-free context info */
439566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
44c4762a1bSJed Brown   m = mctx->m;
45c4762a1bSJed Brown   n = mctx->n;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
489566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
499566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
509566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX1));
519566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX1));
529566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX1));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0, &x0));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX1, &x1));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* dF/dx part */
589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &action));
599566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event1, 0, 0, 0, 0));
60c4762a1bSJed Brown   fos_forward(mctx->tag1, m, n, 0, x0, x1, NULL, action);
61c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
62c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
63c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
6448a46eb9SPierre 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));
65c4762a1bSJed Brown         k++;
66c4762a1bSJed Brown       }
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event1, 0, 0, 0, 0));
70c4762a1bSJed Brown   k = 0;
719566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
729566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));   /*       to INSERT_VALUES and ADD_VALUES         */
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* a * dF/d(xdot) part */
759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event2, 0, 0, 0, 0));
76c4762a1bSJed Brown   fos_forward(mctx->tag2, m, n, 0, x0, x1, NULL, action);
77c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
78c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
79c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
80c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) {
81c4762a1bSJed Brown           action[k] *= mctx->shift;
829566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(Y, 1, &k, &action[k], ADD_VALUES));
83c4762a1bSJed Brown         }
84c4762a1bSJed Brown         k++;
85c4762a1bSJed Brown       }
86c4762a1bSJed Brown     }
87c4762a1bSJed Brown   }
889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event2, 0, 0, 0, 0));
899566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
909566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
919566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Restore local vector */
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX1, &x1));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0, &x0));
969566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX1));
97*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*
101c4762a1bSJed Brown   ADOL-C implementation for Jacobian vector product, using the forward mode of AD.
102c4762a1bSJed Brown   Intended to overload MatMult in matrix-free methods where implicit timestepping
103c4762a1bSJed Brown   has been applied to a problem of the form
104c4762a1bSJed Brown                              du/dt = F(u).
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   Input parameters:
107c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
108c4762a1bSJed Brown   X       - vector to be multiplied by A_shell
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   Output parameters:
111c4762a1bSJed Brown   Y       - product of A_shell and X
112c4762a1bSJed Brown */
113d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianVectorProductIDMass(Mat A_shell, Vec X, Vec Y)
114d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown   AdolcMatCtx       *mctx;
116c4762a1bSJed Brown   PetscInt           m, n, i, j, k = 0, d;
117c4762a1bSJed Brown   const PetscScalar *x0;
118c4762a1bSJed Brown   PetscScalar       *action, *x1;
119c4762a1bSJed Brown   Vec                localX1;
120c4762a1bSJed Brown   DM                 da;
121c4762a1bSJed Brown   DMDALocalInfo      info;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBegin;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Get matrix-free context info */
1269566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
127c4762a1bSJed Brown   m = mctx->m;
128c4762a1bSJed Brown   n = mctx->n;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
1319566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
1329566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
1339566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX1));
1349566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX1));
1359566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX1));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0, &x0));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX1, &x1));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* dF/dx part */
1419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &action));
1429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event1, 0, 0, 0, 0));
143c4762a1bSJed Brown   fos_forward(mctx->tag1, m, n, 0, x0, x1, NULL, action);
144c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
145c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
146c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
14748a46eb9SPierre 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));
148c4762a1bSJed Brown         k++;
149c4762a1bSJed Brown       }
150c4762a1bSJed Brown     }
151c4762a1bSJed Brown   }
1529566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event1, 0, 0, 0, 0));
153c4762a1bSJed Brown   k = 0;
1549566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
1559566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));   /*       to INSERT_VALUES and ADD_VALUES         */
1569566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Restore local vector */
1599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX1, &x1));
1609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0, &x0));
1619566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX1));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* a * dF/d(xdot) part */
1649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event2, 0, 0, 0, 0));
1659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y, mctx->shift, X));
1669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event2, 0, 0, 0, 0));
167*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown   ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
172c4762a1bSJed Brown   Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
173c4762a1bSJed Brown   has been used.
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   Input parameters:
176c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
177c4762a1bSJed Brown   Y       - vector to be multiplied by A_shell transpose
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   Output parameters:
180c4762a1bSJed Brown   X       - product of A_shell transpose and X
181c4762a1bSJed Brown */
182d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianTransposeVectorProduct(Mat A_shell, Vec Y, Vec X)
183d71ae5a4SJacob Faibussowitsch {
184c4762a1bSJed Brown   AdolcMatCtx       *mctx;
185c4762a1bSJed Brown   PetscInt           m, n, i, j, k = 0, d;
186c4762a1bSJed Brown   const PetscScalar *x;
187c4762a1bSJed Brown   PetscScalar       *action, *y;
188c4762a1bSJed Brown   Vec                localY;
189c4762a1bSJed Brown   DM                 da;
190c4762a1bSJed Brown   DMDALocalInfo      info;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   PetscFunctionBegin;
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* Get matrix-free context info */
1959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
196c4762a1bSJed Brown   m = mctx->m;
197c4762a1bSJed Brown   n = mctx->n;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
2009566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
2019566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
2029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localY));
2039566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Y, INSERT_VALUES, localY));
2049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Y, INSERT_VALUES, localY));
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0, &x));
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localY, &y));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* dF/dx part */
2099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &action));
2109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event3, 0, 0, 0, 0));
2119371c9d4SSatish Balay   if (!mctx->flg) zos_forward(mctx->tag1, m, n, 1, x, NULL);
212c4762a1bSJed Brown   fos_reverse(mctx->tag1, m, n, y, action);
213c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
214c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
215c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
21648a46eb9SPierre 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));
217c4762a1bSJed Brown         k++;
218c4762a1bSJed Brown       }
219c4762a1bSJed Brown     }
220c4762a1bSJed Brown   }
2219566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event3, 0, 0, 0, 0));
222c4762a1bSJed Brown   k = 0;
2239566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
2249566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));   /*       to INSERT_VALUES and ADD_VALUES         */
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* a * dF/d(xdot) part */
2279566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event4, 0, 0, 0, 0));
228c4762a1bSJed Brown   if (!mctx->flg) {
229c4762a1bSJed Brown     zos_forward(mctx->tag2, m, n, 1, x, NULL);
230c4762a1bSJed Brown     mctx->flg = PETSC_TRUE;
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown   fos_reverse(mctx->tag2, m, n, y, action);
233c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
234c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
235c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
236c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs + info.xm) && (j >= info.ys) && (j < info.ys + info.ym)) {
237c4762a1bSJed Brown           action[k] *= mctx->shift;
2389566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(X, 1, &k, &action[k], ADD_VALUES));
239c4762a1bSJed Brown         }
240c4762a1bSJed Brown         k++;
241c4762a1bSJed Brown       }
242c4762a1bSJed Brown     }
243c4762a1bSJed Brown   }
2449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event4, 0, 0, 0, 0));
2459566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
2469566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
2479566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Restore local vector */
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localY, &y));
2519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0, &x));
2529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localY));
253*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown /*
257c4762a1bSJed Brown   ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
258c4762a1bSJed Brown   Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
259c4762a1bSJed Brown   has been applied to a problem of the form
260c4762a1bSJed Brown                           du/dt = F(u).
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   Input parameters:
263c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
264c4762a1bSJed Brown   Y       - vector to be multiplied by A_shell transpose
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   Output parameters:
267c4762a1bSJed Brown   X       - product of A_shell transpose and X
268c4762a1bSJed Brown */
269d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscAdolcIJacobianTransposeVectorProductIDMass(Mat A_shell, Vec Y, Vec X)
270d71ae5a4SJacob Faibussowitsch {
271c4762a1bSJed Brown   AdolcMatCtx       *mctx;
272c4762a1bSJed Brown   PetscInt           m, n, i, j, k = 0, d;
273c4762a1bSJed Brown   const PetscScalar *x;
274c4762a1bSJed Brown   PetscScalar       *action, *y;
275c4762a1bSJed Brown   Vec                localY;
276c4762a1bSJed Brown   DM                 da;
277c4762a1bSJed Brown   DMDALocalInfo      info;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   PetscFunctionBegin;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /* Get matrix-free context info */
2829566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
283c4762a1bSJed Brown   m = mctx->m;
284c4762a1bSJed Brown   n = mctx->n;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
2879566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
2889566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
2899566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localY));
2909566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, Y, INSERT_VALUES, localY));
2919566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, Y, INSERT_VALUES, localY));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0, &x));
2939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localY, &y));
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /* dF/dx part */
2969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &action));
2979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event3, 0, 0, 0, 0));
298c4762a1bSJed Brown   if (!mctx->flg) zos_forward(mctx->tag1, m, n, 1, x, NULL);
299c4762a1bSJed Brown   fos_reverse(mctx->tag1, m, n, y, action);
300c4762a1bSJed Brown   for (j = info.gys; j < info.gys + info.gym; j++) {
301c4762a1bSJed Brown     for (i = info.gxs; i < info.gxs + info.gxm; i++) {
302c4762a1bSJed Brown       for (d = 0; d < 2; d++) {
30348a46eb9SPierre 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));
304c4762a1bSJed Brown         k++;
305c4762a1bSJed Brown       }
306c4762a1bSJed Brown     }
307c4762a1bSJed Brown   }
3089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event3, 0, 0, 0, 0));
309c4762a1bSJed Brown   k = 0;
3109566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
3119566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));   /*       to INSERT_VALUES and ADD_VALUES         */
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /* Restore local vector */
3159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localY, &y));
3169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0, &x));
3179566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localY));
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /* a * dF/d(xdot) part */
3209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event4, 0, 0, 0, 0));
3219566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, mctx->shift, Y));
3229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event4, 0, 0, 0, 0));
323*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
324c4762a1bSJed Brown }
325