xref: /petsc/src/ts/tutorials/autodiff/adolc-utils/matfree.cxx (revision 9566063d113dddea24716c546802770db7481bc0)
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 */
30c4762a1bSJed Brown PetscErrorCode PetscAdolcIJacobianVectorProduct(Mat A_shell,Vec X,Vec Y)
31c4762a1bSJed Brown {
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 */
43*9566063dSJacob 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*/
48*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts,&da));
49*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
50*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX1));
51*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX1));
52*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX1));
53c4762a1bSJed Brown 
54*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0,&x0));
55*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX1,&x1));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* dF/dx part */
58*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&action));
59*9566063dSJacob 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++) {
64c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
65*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(Y,1,&k,&action[k],INSERT_VALUES));
66c4762a1bSJed Brown         }
67c4762a1bSJed Brown         k++;
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown     }
70c4762a1bSJed Brown   }
71*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event1,0,0,0,0));
72c4762a1bSJed Brown   k = 0;
73*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
74*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));   /*       to INSERT_VALUES and ADD_VALUES         */
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* a * dF/d(xdot) part */
77*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event2,0,0,0,0));
78c4762a1bSJed Brown   fos_forward(mctx->tag2,m,n,0,x0,x1,NULL,action);
79c4762a1bSJed Brown   for (j=info.gys; j<info.gys+info.gym; j++) {
80c4762a1bSJed Brown     for (i=info.gxs; i<info.gxs+info.gxm; i++) {
81c4762a1bSJed Brown       for (d=0; d<2; d++) {
82c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
83c4762a1bSJed Brown           action[k] *= mctx->shift;
84*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(Y,1,&k,&action[k],ADD_VALUES));
85c4762a1bSJed Brown         }
86c4762a1bSJed Brown         k++;
87c4762a1bSJed Brown       }
88c4762a1bSJed Brown     }
89c4762a1bSJed Brown   }
90*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event2,0,0,0,0));
91*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y));
92*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));
93*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Restore local vector */
96*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX1,&x1));
97*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0,&x0));
98*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX1));
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /*
103c4762a1bSJed Brown   ADOL-C implementation for Jacobian vector product, using the forward mode of AD.
104c4762a1bSJed Brown   Intended to overload MatMult in matrix-free methods where implicit timestepping
105c4762a1bSJed Brown   has been applied to a problem of the form
106c4762a1bSJed Brown                              du/dt = F(u).
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   Input parameters:
109c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
110c4762a1bSJed Brown   X       - vector to be multiplied by A_shell
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   Output parameters:
113c4762a1bSJed Brown   Y       - product of A_shell and X
114c4762a1bSJed Brown */
115c4762a1bSJed Brown PetscErrorCode PetscAdolcIJacobianVectorProductIDMass(Mat A_shell,Vec X,Vec Y)
116c4762a1bSJed Brown {
117c4762a1bSJed Brown   AdolcMatCtx       *mctx;
118c4762a1bSJed Brown   PetscInt          m,n,i,j,k = 0,d;
119c4762a1bSJed Brown   const PetscScalar *x0;
120c4762a1bSJed Brown   PetscScalar       *action,*x1;
121c4762a1bSJed Brown   Vec               localX1;
122c4762a1bSJed Brown   DM                da;
123c4762a1bSJed Brown   DMDALocalInfo     info;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   PetscFunctionBegin;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* Get matrix-free context info */
128*9566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell,&mctx));
129c4762a1bSJed Brown   m = mctx->m;
130c4762a1bSJed Brown   n = mctx->n;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
133*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts,&da));
134*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
135*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX1));
136*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX1));
137*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX1));
138c4762a1bSJed Brown 
139*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0,&x0));
140*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX1,&x1));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* dF/dx part */
143*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&action));
144*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event1,0,0,0,0));
145c4762a1bSJed Brown   fos_forward(mctx->tag1,m,n,0,x0,x1,NULL,action);
146c4762a1bSJed Brown   for (j=info.gys; j<info.gys+info.gym; j++) {
147c4762a1bSJed Brown     for (i=info.gxs; i<info.gxs+info.gxm; i++) {
148c4762a1bSJed Brown       for (d=0; d<2; d++) {
149c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
150*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(Y,1,&k,&action[k],INSERT_VALUES));
151c4762a1bSJed Brown         }
152c4762a1bSJed Brown         k++;
153c4762a1bSJed Brown       }
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown   }
156*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event1,0,0,0,0));
157c4762a1bSJed Brown   k = 0;
158*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Y)); /* Note: Need to assemble between separate calls */
159*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Y));   /*       to INSERT_VALUES and ADD_VALUES         */
160*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* Restore local vector */
163*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX1,&x1));
164*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0,&x0));
165*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX1));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* a * dF/d(xdot) part */
168*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event2,0,0,0,0));
169*9566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y,mctx->shift,X));
170*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event2,0,0,0,0));
171c4762a1bSJed Brown   PetscFunctionReturn(0);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /*
175c4762a1bSJed Brown   ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
176c4762a1bSJed Brown   Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
177c4762a1bSJed Brown   has been used.
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   Input parameters:
180c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
181c4762a1bSJed Brown   Y       - vector to be multiplied by A_shell transpose
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   Output parameters:
184c4762a1bSJed Brown   X       - product of A_shell transpose and X
185c4762a1bSJed Brown */
186c4762a1bSJed Brown PetscErrorCode PetscAdolcIJacobianTransposeVectorProduct(Mat A_shell,Vec Y,Vec X)
187c4762a1bSJed Brown {
188c4762a1bSJed Brown   AdolcMatCtx       *mctx;
189c4762a1bSJed Brown   PetscInt          m,n,i,j,k = 0,d;
190c4762a1bSJed Brown   const PetscScalar *x;
191c4762a1bSJed Brown   PetscScalar       *action,*y;
192c4762a1bSJed Brown   Vec               localY;
193c4762a1bSJed Brown   DM                da;
194c4762a1bSJed Brown   DMDALocalInfo     info;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBegin;
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Get matrix-free context info */
199*9566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell,&mctx));
200c4762a1bSJed Brown   m = mctx->m;
201c4762a1bSJed Brown   n = mctx->n;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
204*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts,&da));
205*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
206*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localY));
207*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,Y,INSERT_VALUES,localY));
208*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,Y,INSERT_VALUES,localY));
209*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0,&x));
210*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localY,&y));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* dF/dx part */
213*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&action));
214*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event3,0,0,0,0));
215c4762a1bSJed Brown   if (!mctx->flg)
216c4762a1bSJed Brown     zos_forward(mctx->tag1,m,n,1,x,NULL);
217c4762a1bSJed Brown   fos_reverse(mctx->tag1,m,n,y,action);
218c4762a1bSJed Brown   for (j=info.gys; j<info.gys+info.gym; j++) {
219c4762a1bSJed Brown     for (i=info.gxs; i<info.gxs+info.gxm; i++) {
220c4762a1bSJed Brown       for (d=0; d<2; d++) {
221c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
222*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(X,1,&k,&action[k],INSERT_VALUES));
223c4762a1bSJed Brown         }
224c4762a1bSJed Brown         k++;
225c4762a1bSJed Brown       }
226c4762a1bSJed Brown     }
227c4762a1bSJed Brown   }
228*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event3,0,0,0,0));
229c4762a1bSJed Brown   k = 0;
230*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
231*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));   /*       to INSERT_VALUES and ADD_VALUES         */
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* a * dF/d(xdot) part */
234*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event4,0,0,0,0));
235c4762a1bSJed Brown   if (!mctx->flg) {
236c4762a1bSJed Brown     zos_forward(mctx->tag2,m,n,1,x,NULL);
237c4762a1bSJed Brown     mctx->flg = PETSC_TRUE;
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown   fos_reverse(mctx->tag2,m,n,y,action);
240c4762a1bSJed Brown   for (j=info.gys; j<info.gys+info.gym; j++) {
241c4762a1bSJed Brown     for (i=info.gxs; i<info.gxs+info.gxm; i++) {
242c4762a1bSJed Brown       for (d=0; d<2; d++) {
243c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
244c4762a1bSJed Brown           action[k] *= mctx->shift;
245*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(X,1,&k,&action[k],ADD_VALUES));
246c4762a1bSJed Brown         }
247c4762a1bSJed Brown         k++;
248c4762a1bSJed Brown       }
249c4762a1bSJed Brown     }
250c4762a1bSJed Brown   }
251*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event4,0,0,0,0));
252*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X));
253*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));
254*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Restore local vector */
257*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localY,&y));
258*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0,&x));
259*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localY));
260c4762a1bSJed Brown   PetscFunctionReturn(0);
261c4762a1bSJed Brown }
262c4762a1bSJed Brown 
263c4762a1bSJed Brown /*
264c4762a1bSJed Brown   ADOL-C implementation for Jacobian transpose vector product, using the reverse mode of AD.
265c4762a1bSJed Brown   Intended to overload MatMultTranspose in matrix-free methods where implicit timestepping
266c4762a1bSJed Brown   has been applied to a problem of the form
267c4762a1bSJed Brown                           du/dt = F(u).
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   Input parameters:
270c4762a1bSJed Brown   A_shell - Jacobian matrix of MatShell type
271c4762a1bSJed Brown   Y       - vector to be multiplied by A_shell transpose
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   Output parameters:
274c4762a1bSJed Brown   X       - product of A_shell transpose and X
275c4762a1bSJed Brown */
276c4762a1bSJed Brown PetscErrorCode PetscAdolcIJacobianTransposeVectorProductIDMass(Mat A_shell,Vec Y,Vec X)
277c4762a1bSJed Brown {
278c4762a1bSJed Brown   AdolcMatCtx       *mctx;
279c4762a1bSJed Brown   PetscInt          m,n,i,j,k = 0,d;
280c4762a1bSJed Brown   const PetscScalar *x;
281c4762a1bSJed Brown   PetscScalar       *action,*y;
282c4762a1bSJed Brown   Vec               localY;
283c4762a1bSJed Brown   DM                da;
284c4762a1bSJed Brown   DMDALocalInfo     info;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   PetscFunctionBegin;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   /* Get matrix-free context info */
289*9566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell,&mctx));
290c4762a1bSJed Brown   m = mctx->m;
291c4762a1bSJed Brown   n = mctx->n;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   /* Get local input vectors and extract data, x0 and x1*/
294*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts,&da));
295*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
296*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localY));
297*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,Y,INSERT_VALUES,localY));
298*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,Y,INSERT_VALUES,localY));
299*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mctx->localX0,&x));
300*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localY,&y));
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /* dF/dx part */
303*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&action));
304*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event3,0,0,0,0));
305c4762a1bSJed Brown   if (!mctx->flg) zos_forward(mctx->tag1,m,n,1,x,NULL);
306c4762a1bSJed Brown   fos_reverse(mctx->tag1,m,n,y,action);
307c4762a1bSJed Brown   for (j=info.gys; j<info.gys+info.gym; j++) {
308c4762a1bSJed Brown     for (i=info.gxs; i<info.gxs+info.gxm; i++) {
309c4762a1bSJed Brown       for (d=0; d<2; d++) {
310c4762a1bSJed Brown         if ((i >= info.xs) && (i < info.xs+info.xm) && (j >= info.ys) && (j < info.ys+info.ym)) {
311*9566063dSJacob Faibussowitsch           PetscCall(VecSetValuesLocal(X,1,&k,&action[k],INSERT_VALUES));
312c4762a1bSJed Brown         }
313c4762a1bSJed Brown         k++;
314c4762a1bSJed Brown       }
315c4762a1bSJed Brown     }
316c4762a1bSJed Brown   }
317*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event3,0,0,0,0));
318c4762a1bSJed Brown   k = 0;
319*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(X)); /* Note: Need to assemble between separate calls */
320*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(X));   /*       to INSERT_VALUES and ADD_VALUES         */
321*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(action));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* Restore local vector */
324*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localY,&y));
325*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mctx->localX0,&x));
326*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localY));
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /* a * dF/d(xdot) part */
329*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(mctx->event4,0,0,0,0));
330*9566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,mctx->shift,Y));
331*9566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(mctx->event4,0,0,0,0));
332c4762a1bSJed Brown   PetscFunctionReturn(0);
333c4762a1bSJed Brown }
334