xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5opt_ic.c (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    See ex5.c for details on the equation.
5c4762a1bSJed Brown    This code demonestrates the TSAdjoint/TAO interface to solve an inverse initial value problem built on a system of
6c4762a1bSJed Brown    time-dependent partial differential equations.
7c4762a1bSJed Brown    In this problem, the initial value for the PDE is unknown, but the output (the final solution of the PDE) is known.
8c4762a1bSJed Brown    We want to determine the initial value that can produce the given output.
9c4762a1bSJed Brown    We formulate the problem as a nonlinear optimization problem that minimizes the discrepency between the simulated
10c4762a1bSJed Brown    result and given reference solution, calculate the gradient of the objective function with the discrete adjoint
11c4762a1bSJed Brown    solver, and solve the optimization problem with TAO.
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    Runtime options:
14c4762a1bSJed Brown      -forwardonly  - run only the forward simulation
15c4762a1bSJed Brown      -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used
16c4762a1bSJed Brown  */
17c4762a1bSJed Brown 
18c4762a1bSJed Brown #include <petscsys.h>
19c4762a1bSJed Brown #include <petscdm.h>
20c4762a1bSJed Brown #include <petscdmda.h>
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown #include <petsctao.h>
23c4762a1bSJed Brown 
24c4762a1bSJed Brown typedef struct {
25c4762a1bSJed Brown   PetscScalar u,v;
26c4762a1bSJed Brown } Field;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   PetscReal D1,D2,gamma,kappa;
30c4762a1bSJed Brown   TS        ts;
31c4762a1bSJed Brown   Vec       U;
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /* User-defined routines */
35c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*),InitialConditions(DM,Vec);
36c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
37c4762a1bSJed Brown extern PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
38c4762a1bSJed Brown extern PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
39c4762a1bSJed Brown extern PetscErrorCode FormFunctionAndGradient(Tao,Vec,PetscReal*,Vec,void*);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*
42c4762a1bSJed Brown    Set terminal condition for the adjoint variable
43c4762a1bSJed Brown  */
44c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,Vec U,AppCtx *appctx)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN]="";
47c4762a1bSJed Brown   PetscViewer    viewer;
48c4762a1bSJed Brown   Vec            Uob;
49c4762a1bSJed Brown   PetscErrorCode ierr;
50c4762a1bSJed Brown 
51*362febeeSStefano Zampini   PetscFunctionBegin;
52c4762a1bSJed Brown   ierr = VecDuplicate(U,&Uob);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = VecLoad(Uob,viewer);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = VecAYPX(Uob,-1.,U);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = VecScale(Uob,2.0);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = VecAXPY(lambda,1.,Uob);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = VecDestroy(&Uob);CHKERRQ(ierr);
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*
65c4762a1bSJed Brown    Set up a viewer that dumps data into a binary file
66c4762a1bSJed Brown  */
67c4762a1bSJed Brown PetscErrorCode OutputBIN(DM da, const char *filename, PetscViewer *viewer)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscErrorCode ierr;
70c4762a1bSJed Brown 
71*362febeeSStefano Zampini   PetscFunctionBegin;
72c4762a1bSJed Brown   ierr = PetscViewerCreate(PetscObjectComm((PetscObject)da),viewer);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
76c4762a1bSJed Brown   PetscFunctionReturn(0);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown /*
80c4762a1bSJed Brown    Generate a reference solution and save it to a binary file
81c4762a1bSJed Brown  */
82c4762a1bSJed Brown PetscErrorCode GenerateOBs(TS ts,Vec U,AppCtx *appctx)
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN] = "";
85c4762a1bSJed Brown   PetscViewer    viewer;
86c4762a1bSJed Brown   DM             da;
87c4762a1bSJed Brown   PetscErrorCode ierr;
88c4762a1bSJed Brown 
89*362febeeSStefano Zampini   PetscFunctionBegin;
90c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = OutputBIN(da,filename,&viewer);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = VecView(U,viewer);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
100c4762a1bSJed Brown {
101c4762a1bSJed Brown   PetscErrorCode ierr;
102c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
103c4762a1bSJed Brown   Field          **u;
104c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   PetscFunctionBegin;
107c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
110c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
113c4762a1bSJed Brown   /* Get local grid boundaries */
114c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
117c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
118c4762a1bSJed Brown     y = j*hy;
119c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
120c4762a1bSJed Brown       x = i*hx;
121c4762a1bSJed Brown       if ((1.0 <= x) && (x <= 1.5) && (1.0 <= y) && (y <= 1.5)) u[j][i].v = .25*PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0);
122c4762a1bSJed Brown       else u[j][i].v = 0.0;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
125c4762a1bSJed Brown     }
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
129c4762a1bSJed Brown   PetscFunctionReturn(0);
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions(DM da,Vec U)
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   PetscErrorCode ierr;
135c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
136c4762a1bSJed Brown   Field          **u;
137c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   PetscFunctionBegin;
140c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
143c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
146c4762a1bSJed Brown   /* Get local grid boundaries */
147c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
150c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
151c4762a1bSJed Brown     y = j*hy;
152c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
153c4762a1bSJed Brown       x = i*hx;
154c4762a1bSJed Brown       if ((1.0 <= x) && (x <= 1.5) && (1.0 <= y) && (y <= 1.5)) u[j][i].v = .25*PetscPowReal(PetscSinReal(4.0*PETSC_PI*0.5*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*0.5*y),2.0);
155c4762a1bSJed Brown       else u[j][i].v = 0.0;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
158c4762a1bSJed Brown     }
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
162c4762a1bSJed Brown   PetscFunctionReturn(0);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions2(DM da,Vec U)
166c4762a1bSJed Brown {
167c4762a1bSJed Brown   PetscErrorCode ierr;
168c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
169c4762a1bSJed Brown   Field          **u;
170c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBegin;
173c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
176c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
179c4762a1bSJed Brown   /* Get local grid boundaries */
180c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
183c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
184c4762a1bSJed Brown     y = j*hy;
185c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
186c4762a1bSJed Brown       x = i*hx;
187c803a25aSBarry Smith       if (PetscGTE(x,1.0) && PetscLTE(x,1.5) && PetscGTE(y,1.0) && PetscLTE(y,1.5)) u[j][i].v = PetscPowReal(PetscSinReal(4.0*PETSC_PI*(x-0.5)),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0)/4.0;
188c4762a1bSJed Brown       else u[j][i].v = 0.0;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
191c4762a1bSJed Brown     }
192c4762a1bSJed Brown   }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
195c4762a1bSJed Brown   PetscFunctionReturn(0);
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions3(DM da,Vec U)
199c4762a1bSJed Brown {
200c4762a1bSJed Brown   PetscErrorCode ierr;
201c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
202c4762a1bSJed Brown   Field          **u;
203c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBegin;
206c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
209c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
212c4762a1bSJed Brown   /* Get local grid boundaries */
213c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
216c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
217c4762a1bSJed Brown     y = j*hy;
218c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
219c4762a1bSJed Brown       x = i*hx;
220c4762a1bSJed Brown       if ((0.5 <= x) && (x <= 2.0) && (0.5 <= y) && (y <= 2.0)) u[j][i].v = .25*PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0);
221c4762a1bSJed Brown       else u[j][i].v = 0.0;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
224c4762a1bSJed Brown     }
225c4762a1bSJed Brown   }
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
228c4762a1bSJed Brown   PetscFunctionReturn(0);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown int main(int argc,char **argv)
232c4762a1bSJed Brown {
233c4762a1bSJed Brown   PetscErrorCode ierr;
234c4762a1bSJed Brown   DM             da;
235c4762a1bSJed Brown   AppCtx         appctx;
236c4762a1bSJed Brown   PetscBool      forwardonly = PETSC_FALSE,implicitform = PETSC_FALSE;
237c4762a1bSJed Brown   PetscInt       perturbic = 1;
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
240c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr);
241c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-perturbic",&perturbic,NULL);CHKERRQ(ierr);
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
245c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
246c4762a1bSJed Brown   appctx.gamma = .024;
247c4762a1bSJed Brown   appctx.kappa = .06;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Create distributed array (DMDA) to manage parallel grid and vectors */
250c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
253c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
254c4762a1bSJed Brown   ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr);
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Extract global vectors from DMDA; then duplicate for remaining
257c4762a1bSJed Brown      vectors that are the same types */
258c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&appctx.U);CHKERRQ(ierr);
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* Create timestepping solver context */
261c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&appctx.ts);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = TSSetType(appctx.ts,TSCN);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = TSSetDM(appctx.ts,da);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = TSSetProblemType(appctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = TSSetEquationType(appctx.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
266c4762a1bSJed Brown   if (!implicitform) {
267c4762a1bSJed Brown     ierr = TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr);
268c4762a1bSJed Brown     ierr = TSSetRHSJacobian(appctx.ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr);
269c4762a1bSJed Brown   } else {
270c4762a1bSJed Brown     ierr = TSSetIFunction(appctx.ts,NULL,IFunction,&appctx);CHKERRQ(ierr);
271c4762a1bSJed Brown     ierr = TSSetIJacobian(appctx.ts,NULL,NULL,IJacobian,&appctx);CHKERRQ(ierr);
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* Set initial conditions */
275c4762a1bSJed Brown   ierr = InitialConditions(da,appctx.U);CHKERRQ(ierr);
276c4762a1bSJed Brown   ierr = TSSetSolution(appctx.ts,appctx.U);CHKERRQ(ierr);
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /* Set solver options */
279c4762a1bSJed Brown   ierr = TSSetMaxTime(appctx.ts,2000.0);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = TSSetTimeStep(appctx.ts,0.5);CHKERRQ(ierr);
281c4762a1bSJed Brown   ierr = TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
282c4762a1bSJed Brown   ierr = TSSetFromOptions(appctx.ts);CHKERRQ(ierr);
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   ierr = GenerateOBs(appctx.ts,appctx.U,&appctx);CHKERRQ(ierr);
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   if (!forwardonly) {
287c4762a1bSJed Brown     Tao           tao;
288c4762a1bSJed Brown     Vec           P;
289c4762a1bSJed Brown     Vec           lambda[1];
290956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
291c4762a1bSJed Brown     PetscLogStage opt_stage;
292956f8c0dSBarry Smith #endif
293c4762a1bSJed Brown 
294c4762a1bSJed Brown     ierr = PetscLogStageRegister("Optimization",&opt_stage);CHKERRQ(ierr);
295c4762a1bSJed Brown     ierr = PetscLogStagePush(opt_stage);CHKERRQ(ierr);
296c4762a1bSJed Brown     if (perturbic == 1) {
297c4762a1bSJed Brown       ierr = PerturbedInitialConditions(da,appctx.U);CHKERRQ(ierr);
298c4762a1bSJed Brown     } else if (perturbic == 2) {
299c4762a1bSJed Brown       ierr = PerturbedInitialConditions2(da,appctx.U);CHKERRQ(ierr);
300c4762a1bSJed Brown     } else if (perturbic == 3) {
301c4762a1bSJed Brown       ierr = PerturbedInitialConditions3(da,appctx.U);CHKERRQ(ierr);
302c4762a1bSJed Brown     }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     ierr = VecDuplicate(appctx.U,&lambda[0]);CHKERRQ(ierr);
305c4762a1bSJed Brown     ierr = TSSetCostGradients(appctx.ts,1,lambda,NULL);CHKERRQ(ierr);
306c4762a1bSJed Brown 
307c4762a1bSJed Brown     /* Have the TS save its trajectory needed by TSAdjointSolve() */
308c4762a1bSJed Brown     ierr = TSSetSaveTrajectory(appctx.ts);CHKERRQ(ierr);
309c4762a1bSJed Brown 
310c4762a1bSJed Brown     /* Create TAO solver and set desired solution method */
311c4762a1bSJed Brown     ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
312c4762a1bSJed Brown     ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
313c4762a1bSJed Brown 
314c4762a1bSJed Brown     /* Set initial guess for TAO */
315c4762a1bSJed Brown     ierr = VecDuplicate(appctx.U,&P);CHKERRQ(ierr);
316c4762a1bSJed Brown     ierr = VecCopy(appctx.U,P);CHKERRQ(ierr);
317c4762a1bSJed Brown     ierr = TaoSetInitialVector(tao,P);CHKERRQ(ierr);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown     /* Set routine for function and gradient evaluation */
320c4762a1bSJed Brown     ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionAndGradient,&appctx);CHKERRQ(ierr);
321c4762a1bSJed Brown 
322c4762a1bSJed Brown     /* Check for any TAO command line options */
323c4762a1bSJed Brown     ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     ierr = TaoSolve(tao);CHKERRQ(ierr);
326c4762a1bSJed Brown     ierr = TaoDestroy(&tao);CHKERRQ(ierr);
327c4762a1bSJed Brown     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
328c4762a1bSJed Brown     ierr = VecDestroy(&P);CHKERRQ(ierr);
329c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
330c4762a1bSJed Brown   }
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /* Free work space.  All PETSc objects should be destroyed when they
333c4762a1bSJed Brown      are no longer needed. */
334c4762a1bSJed Brown   ierr = VecDestroy(&appctx.U);CHKERRQ(ierr);
335c4762a1bSJed Brown   ierr = TSDestroy(&appctx.ts);CHKERRQ(ierr);
336c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
337c4762a1bSJed Brown   ierr = PetscFinalize();
338c4762a1bSJed Brown   return ierr;
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown /* ------------------------ TS callbacks ---------------------------- */
342c4762a1bSJed Brown 
343c4762a1bSJed Brown /*
344c4762a1bSJed Brown    RHSFunction - Evaluates nonlinear function, F(x).
345c4762a1bSJed Brown 
346c4762a1bSJed Brown    Input Parameters:
347c4762a1bSJed Brown .  ts - the TS context
348c4762a1bSJed Brown .  X - input vector
349c4762a1bSJed Brown .  ptr - optional user-defined context, as set by TSSetRHSFunction()
350c4762a1bSJed Brown 
351c4762a1bSJed Brown    Output Parameter:
352c4762a1bSJed Brown .  F - function vector
353c4762a1bSJed Brown  */
354c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
355c4762a1bSJed Brown {
356c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
357c4762a1bSJed Brown   DM             da;
358c4762a1bSJed Brown   PetscErrorCode ierr;
359c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
360c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
361c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
362c4762a1bSJed Brown   Field          **u,**f;
363c4762a1bSJed Brown   Vec            localU;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   PetscFunctionBegin;
366c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
369c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
370c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* Scatter ghost points to local vector,using the 2-step process
373c4762a1bSJed Brown      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
374c4762a1bSJed Brown      By placing code between these two statements, computations can be
375c4762a1bSJed Brown      done while messages are in transition. */
376c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
378c4762a1bSJed Brown 
379c4762a1bSJed Brown   /* Get pointers to vector data */
380c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
381c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   /* Get local grid boundaries */
384c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
387c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
388c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
389c4762a1bSJed Brown       uc        = u[j][i].u;
390c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
391c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
392c4762a1bSJed Brown       vc        = u[j][i].v;
393c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
394c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
395c4762a1bSJed Brown       f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
396c4762a1bSJed Brown       f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
397c4762a1bSJed Brown     }
398c4762a1bSJed Brown   }
399ca0c957dSBarry Smith   ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr);
400c4762a1bSJed Brown 
401c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
402c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
403c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
404c4762a1bSJed Brown   PetscFunctionReturn(0);
405c4762a1bSJed Brown }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx)
408c4762a1bSJed Brown {
409c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
410c4762a1bSJed Brown   DM             da;
411c4762a1bSJed Brown   PetscErrorCode ierr;
412c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
413c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
414c4762a1bSJed Brown   PetscScalar    uc,vc;
415c4762a1bSJed Brown   Field          **u;
416c4762a1bSJed Brown   Vec            localU;
417c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
418c4762a1bSJed Brown   PetscScalar    entries[6];
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   PetscFunctionBegin;
421c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
422c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
423c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
426c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
427c4762a1bSJed Brown 
428c4762a1bSJed Brown   /* Scatter ghost points to local vector,using the 2-step process
429c4762a1bSJed Brown      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
430c4762a1bSJed Brown      By placing code between these two statements, computations can be
431c4762a1bSJed Brown      done while messages are in transition. */
432c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
433c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   /* Get pointers to vector data */
436c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /* Get local grid boundaries */
439c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   stencil[0].k = 0;
442c4762a1bSJed Brown   stencil[1].k = 0;
443c4762a1bSJed Brown   stencil[2].k = 0;
444c4762a1bSJed Brown   stencil[3].k = 0;
445c4762a1bSJed Brown   stencil[4].k = 0;
446c4762a1bSJed Brown   stencil[5].k = 0;
447c4762a1bSJed Brown   rowstencil.k = 0;
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
450c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
451c4762a1bSJed Brown     stencil[0].j = j-1;
452c4762a1bSJed Brown     stencil[1].j = j+1;
453c4762a1bSJed Brown     stencil[2].j = j;
454c4762a1bSJed Brown     stencil[3].j = j;
455c4762a1bSJed Brown     stencil[4].j = j;
456c4762a1bSJed Brown     stencil[5].j = j;
457c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
458c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
459c4762a1bSJed Brown       uc = u[j][i].u;
460c4762a1bSJed Brown       vc = u[j][i].v;
461c4762a1bSJed Brown 
462c4762a1bSJed Brown       /* uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
463c4762a1bSJed Brown          uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
464c4762a1bSJed Brown          vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
465c4762a1bSJed Brown          vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
466c4762a1bSJed Brown          f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
467c4762a1bSJed Brown 
468c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
469c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
470c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
471c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
472c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
473c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
474c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
475c4762a1bSJed Brown 
476c4762a1bSJed Brown       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
477c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = appctx->D2*sy;
478c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = appctx->D2*sy;
479c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = appctx->D2*sx;
480c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = appctx->D2*sx;
481c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
482c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = vc*vc;
483c4762a1bSJed Brown       rowstencil.c = 1;
484c4762a1bSJed Brown 
485c4762a1bSJed Brown       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
486c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
487c4762a1bSJed Brown     }
488c4762a1bSJed Brown   }
489ca0c957dSBarry Smith   ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr);
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
492c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
493c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
494c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
495c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
496c4762a1bSJed Brown   PetscFunctionReturn(0);
497c4762a1bSJed Brown }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown /*
500c4762a1bSJed Brown    IFunction - Evaluates implicit nonlinear function, xdot - F(x).
501c4762a1bSJed Brown 
502c4762a1bSJed Brown    Input Parameters:
503c4762a1bSJed Brown .  ts - the TS context
504c4762a1bSJed Brown .  U - input vector
505c4762a1bSJed Brown .  Udot - input vector
506c4762a1bSJed Brown .  ptr - optional user-defined context, as set by TSSetRHSFunction()
507c4762a1bSJed Brown 
508c4762a1bSJed Brown    Output Parameter:
509c4762a1bSJed Brown .  F - function vector
510c4762a1bSJed Brown  */
511c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
512c4762a1bSJed Brown {
513c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
514c4762a1bSJed Brown   DM             da;
515c4762a1bSJed Brown   PetscErrorCode ierr;
516c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
517c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
518c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
519c4762a1bSJed Brown   Field          **u,**f,**udot;
520c4762a1bSJed Brown   Vec            localU;
521c4762a1bSJed Brown 
522c4762a1bSJed Brown   PetscFunctionBegin;
523c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
524c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
525c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
526c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
527c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   /* Scatter ghost points to local vector,using the 2-step process
530c4762a1bSJed Brown      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
531c4762a1bSJed Brown      By placing code between these two statements, computations can be
532c4762a1bSJed Brown      done while messages are in transition. */
533c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
534c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
535c4762a1bSJed Brown 
536c4762a1bSJed Brown   /* Get pointers to vector data */
537c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
538c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
539c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr);
540c4762a1bSJed Brown 
541c4762a1bSJed Brown   /* Get local grid boundaries */
542c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
543c4762a1bSJed Brown 
544c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
545c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
546c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
547c4762a1bSJed Brown       uc        = u[j][i].u;
548c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
549c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
550c4762a1bSJed Brown       vc        = u[j][i].v;
551c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
552c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
553c4762a1bSJed Brown       f[j][i].u = udot[j][i].u - ( appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc));
554c4762a1bSJed Brown       f[j][i].v = udot[j][i].v - ( appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc);
555c4762a1bSJed Brown     }
556c4762a1bSJed Brown   }
557ca0c957dSBarry Smith   ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr);
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
560c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
561c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr);
562c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
563c4762a1bSJed Brown   PetscFunctionReturn(0);
564c4762a1bSJed Brown }
565c4762a1bSJed Brown 
566c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx)
567c4762a1bSJed Brown {
568c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
569c4762a1bSJed Brown   DM             da;
570c4762a1bSJed Brown   PetscErrorCode ierr;
571c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
572c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
573c4762a1bSJed Brown   PetscScalar    uc,vc;
574c4762a1bSJed Brown   Field          **u;
575c4762a1bSJed Brown   Vec            localU;
576c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
577c4762a1bSJed Brown   PetscScalar    entries[6];
578c4762a1bSJed Brown 
579c4762a1bSJed Brown   PetscFunctionBegin;
580c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
581c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
582c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
583c4762a1bSJed Brown 
584c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
585c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
586c4762a1bSJed Brown 
587c4762a1bSJed Brown   /* Scatter ghost points to local vector,using the 2-step process
588c4762a1bSJed Brown      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
589c4762a1bSJed Brown      By placing code between these two statements, computations can be
590c4762a1bSJed Brown      done while messages are in transition.*/
591c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
592c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   /* Get pointers to vector data */
595c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
596c4762a1bSJed Brown 
597c4762a1bSJed Brown   /* Get local grid boundaries */
598c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
599c4762a1bSJed Brown 
600c4762a1bSJed Brown   stencil[0].k = 0;
601c4762a1bSJed Brown   stencil[1].k = 0;
602c4762a1bSJed Brown   stencil[2].k = 0;
603c4762a1bSJed Brown   stencil[3].k = 0;
604c4762a1bSJed Brown   stencil[4].k = 0;
605c4762a1bSJed Brown   stencil[5].k = 0;
606c4762a1bSJed Brown   rowstencil.k = 0;
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
609c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
610c4762a1bSJed Brown     stencil[0].j = j-1;
611c4762a1bSJed Brown     stencil[1].j = j+1;
612c4762a1bSJed Brown     stencil[2].j = j;
613c4762a1bSJed Brown     stencil[3].j = j;
614c4762a1bSJed Brown     stencil[4].j = j;
615c4762a1bSJed Brown     stencil[5].j = j;
616c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
617c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
618c4762a1bSJed Brown       uc = u[j][i].u;
619c4762a1bSJed Brown       vc = u[j][i].v;
620c4762a1bSJed Brown 
621c4762a1bSJed Brown       /* uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
622c4762a1bSJed Brown          uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
623c4762a1bSJed Brown          vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
624c4762a1bSJed Brown          vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
625c4762a1bSJed Brown          f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc); */
626c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy;
627c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy;
628c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx;
629c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx;
630c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a;
631c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc;
632c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
633c4762a1bSJed Brown 
634c4762a1bSJed Brown       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
635c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = -appctx->D2*sy;
636c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = -appctx->D2*sy;
637c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = -appctx->D2*sx;
638c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = -appctx->D2*sx;
639c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a;
640c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = -vc*vc;
641c4762a1bSJed Brown       rowstencil.c = 1;
642c4762a1bSJed Brown 
643c4762a1bSJed Brown       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
644c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
645c4762a1bSJed Brown     }
646c4762a1bSJed Brown   }
647ca0c957dSBarry Smith   ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr);
648c4762a1bSJed Brown 
649c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
650c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
651c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
652c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
653c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
654c4762a1bSJed Brown   PetscFunctionReturn(0);
655c4762a1bSJed Brown }
656c4762a1bSJed Brown 
657c4762a1bSJed Brown /* ------------------------ TAO callbacks ---------------------------- */
658c4762a1bSJed Brown 
659c4762a1bSJed Brown /*
660c4762a1bSJed Brown    FormFunctionAndGradient - Evaluates the function and corresponding gradient.
661c4762a1bSJed Brown 
662c4762a1bSJed Brown    Input Parameters:
663c4762a1bSJed Brown    tao - the Tao context
664c4762a1bSJed Brown    P   - the input vector
665c4762a1bSJed Brown    ctx - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
666c4762a1bSJed Brown 
667c4762a1bSJed Brown    Output Parameters:
668c4762a1bSJed Brown    f   - the newly evaluated function
669c4762a1bSJed Brown    G   - the newly evaluated gradient
670c4762a1bSJed Brown */
671c4762a1bSJed Brown PetscErrorCode FormFunctionAndGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
672c4762a1bSJed Brown {
673c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
674c4762a1bSJed Brown   PetscReal      soberr,timestep;
675c4762a1bSJed Brown   Vec            *lambda;
676c4762a1bSJed Brown   Vec            SDiff;
677c4762a1bSJed Brown   DM             da;
678c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN]="";
679c4762a1bSJed Brown   PetscViewer    viewer;
680c4762a1bSJed Brown   PetscErrorCode ierr;
681c4762a1bSJed Brown 
682c4762a1bSJed Brown   PetscFunctionBeginUser;
683c4762a1bSJed Brown   ierr = TSSetTime(appctx->ts,0.0);CHKERRQ(ierr);
684c4762a1bSJed Brown   ierr = TSGetTimeStep(appctx->ts,&timestep);CHKERRQ(ierr);
685c4762a1bSJed Brown   if (timestep<0) {
686c4762a1bSJed Brown     ierr = TSSetTimeStep(appctx->ts,-timestep);CHKERRQ(ierr);
687c4762a1bSJed Brown   }
688c4762a1bSJed Brown   ierr = TSSetStepNumber(appctx->ts,0);CHKERRQ(ierr);
689c4762a1bSJed Brown   ierr = TSSetFromOptions(appctx->ts);CHKERRQ(ierr);
690c4762a1bSJed Brown 
691c4762a1bSJed Brown   ierr = VecDuplicate(P,&SDiff);CHKERRQ(ierr);
692c4762a1bSJed Brown   ierr = VecCopy(P,appctx->U);CHKERRQ(ierr);
693c4762a1bSJed Brown   ierr = TSGetDM(appctx->ts,&da);CHKERRQ(ierr);
694c4762a1bSJed Brown   *f = 0;
695c4762a1bSJed Brown 
696c4762a1bSJed Brown   ierr = TSSolve(appctx->ts,appctx->U);CHKERRQ(ierr);
697c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
698c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
699c4762a1bSJed Brown   ierr = VecLoad(SDiff,viewer);CHKERRQ(ierr);
700c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
701c4762a1bSJed Brown   ierr = VecAYPX(SDiff,-1.,appctx->U);CHKERRQ(ierr);
702c4762a1bSJed Brown   ierr = VecDot(SDiff,SDiff,&soberr);CHKERRQ(ierr);
703c4762a1bSJed Brown   *f += soberr;
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   ierr = TSGetCostGradients(appctx->ts,NULL,&lambda,NULL);CHKERRQ(ierr);
706c4762a1bSJed Brown   ierr = VecSet(lambda[0],0.0);CHKERRQ(ierr);
707c4762a1bSJed Brown   ierr = InitializeLambda(da,lambda[0],appctx->U,appctx);CHKERRQ(ierr);
708c4762a1bSJed Brown   ierr = TSAdjointSolve(appctx->ts);CHKERRQ(ierr);
709c4762a1bSJed Brown 
710c4762a1bSJed Brown   ierr = VecCopy(lambda[0],G);CHKERRQ(ierr);
711c4762a1bSJed Brown 
712c4762a1bSJed Brown   ierr = VecDestroy(&SDiff);CHKERRQ(ierr);
713c4762a1bSJed Brown   PetscFunctionReturn(0);
714c4762a1bSJed Brown }
715c4762a1bSJed Brown 
716c4762a1bSJed Brown /*TEST
717c4762a1bSJed Brown 
718c4762a1bSJed Brown    build:
719c4762a1bSJed Brown       requires: !complex !single
720c4762a1bSJed Brown 
721c4762a1bSJed Brown    test:
722c4762a1bSJed Brown       args: -ts_max_steps 5 -ts_type rk -ts_rk_type 3 -ts_trajectory_type memory -tao_monitor -tao_view -tao_gatol 1e-6
723c4762a1bSJed Brown       output_file: output/ex5opt_ic_1.out
724c4762a1bSJed Brown 
725c4762a1bSJed Brown TEST*/
726