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