xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5opt_ic.c (revision 60f0b76ed2293387b137bb632a681b21c2f7fb7d)
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 
18*60f0b76eSHong Zhang #include "reaction_diffusion.h"
19c4762a1bSJed Brown #include <petscdm.h>
20c4762a1bSJed Brown #include <petscdmda.h>
21c4762a1bSJed Brown #include <petsctao.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /* User-defined routines */
24c4762a1bSJed Brown extern PetscErrorCode FormFunctionAndGradient(Tao,Vec,PetscReal*,Vec,void*);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown    Set terminal condition for the adjoint variable
28c4762a1bSJed Brown  */
29c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,Vec U,AppCtx *appctx)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN]="";
32c4762a1bSJed Brown   PetscViewer    viewer;
33c4762a1bSJed Brown   Vec            Uob;
34c4762a1bSJed Brown   PetscErrorCode ierr;
35c4762a1bSJed Brown 
36362febeeSStefano Zampini   PetscFunctionBegin;
37c4762a1bSJed Brown   ierr = VecDuplicate(U,&Uob);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = VecLoad(Uob,viewer);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecAYPX(Uob,-1.,U);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecScale(Uob,2.0);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = VecAXPY(lambda,1.,Uob);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = VecDestroy(&Uob);CHKERRQ(ierr);
46c4762a1bSJed Brown   PetscFunctionReturn(0);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown    Set up a viewer that dumps data into a binary file
51c4762a1bSJed Brown  */
52c4762a1bSJed Brown PetscErrorCode OutputBIN(DM da, const char *filename, PetscViewer *viewer)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   PetscErrorCode ierr;
55c4762a1bSJed Brown 
56362febeeSStefano Zampini   PetscFunctionBegin;
57c4762a1bSJed Brown   ierr = PetscViewerCreate(PetscObjectComm((PetscObject)da),viewer);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*
65c4762a1bSJed Brown    Generate a reference solution and save it to a binary file
66c4762a1bSJed Brown  */
67c4762a1bSJed Brown PetscErrorCode GenerateOBs(TS ts,Vec U,AppCtx *appctx)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN] = "";
70c4762a1bSJed Brown   PetscViewer    viewer;
71c4762a1bSJed Brown   DM             da;
72c4762a1bSJed Brown   PetscErrorCode ierr;
73c4762a1bSJed Brown 
74362febeeSStefano Zampini   PetscFunctionBegin;
75c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = OutputBIN(da,filename,&viewer);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = VecView(U,viewer);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
85c4762a1bSJed Brown {
86c4762a1bSJed Brown   PetscErrorCode ierr;
87c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
88c4762a1bSJed Brown   Field          **u;
89c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   PetscFunctionBegin;
92c4762a1bSJed 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);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
95c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
98c4762a1bSJed Brown   /* Get local grid boundaries */
99c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
102c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
103c4762a1bSJed Brown     y = j*hy;
104c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
105c4762a1bSJed Brown       x = i*hx;
106c4762a1bSJed 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);
107c4762a1bSJed Brown       else u[j][i].v = 0.0;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
110c4762a1bSJed Brown     }
111c4762a1bSJed Brown   }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
114c4762a1bSJed Brown   PetscFunctionReturn(0);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions(DM da,Vec U)
118c4762a1bSJed Brown {
119c4762a1bSJed Brown   PetscErrorCode ierr;
120c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
121c4762a1bSJed Brown   Field          **u;
122c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBegin;
125c4762a1bSJed 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);
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
128c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
131c4762a1bSJed Brown   /* Get local grid boundaries */
132c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
135c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
136c4762a1bSJed Brown     y = j*hy;
137c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
138c4762a1bSJed Brown       x = i*hx;
139c4762a1bSJed 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);
140c4762a1bSJed Brown       else u[j][i].v = 0.0;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
143c4762a1bSJed Brown     }
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
147c4762a1bSJed Brown   PetscFunctionReturn(0);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions2(DM da,Vec U)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   PetscErrorCode ierr;
153c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
154c4762a1bSJed Brown   Field          **u;
155c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBegin;
158c4762a1bSJed 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);
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
161c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
164c4762a1bSJed Brown   /* Get local grid boundaries */
165c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
168c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
169c4762a1bSJed Brown     y = j*hy;
170c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
171c4762a1bSJed Brown       x = i*hx;
172c803a25aSBarry 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;
173c4762a1bSJed Brown       else u[j][i].v = 0.0;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
176c4762a1bSJed Brown     }
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
180c4762a1bSJed Brown   PetscFunctionReturn(0);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown PetscErrorCode PerturbedInitialConditions3(DM da,Vec U)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   PetscErrorCode ierr;
186c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
187c4762a1bSJed Brown   Field          **u;
188c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBegin;
191c4762a1bSJed 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);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
194c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
197c4762a1bSJed Brown   /* Get local grid boundaries */
198c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
201c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
202c4762a1bSJed Brown     y = j*hy;
203c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
204c4762a1bSJed Brown       x = i*hx;
205c4762a1bSJed 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);
206c4762a1bSJed Brown       else u[j][i].v = 0.0;
207c4762a1bSJed Brown 
208c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
209c4762a1bSJed Brown     }
210c4762a1bSJed Brown   }
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
213c4762a1bSJed Brown   PetscFunctionReturn(0);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown int main(int argc,char **argv)
217c4762a1bSJed Brown {
218c4762a1bSJed Brown   PetscErrorCode ierr;
219c4762a1bSJed Brown   DM             da;
220c4762a1bSJed Brown   AppCtx         appctx;
221c4762a1bSJed Brown   PetscBool      forwardonly = PETSC_FALSE,implicitform = PETSC_FALSE;
222c4762a1bSJed Brown   PetscInt       perturbic = 1;
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
225c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr);
227c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-perturbic",&perturbic,NULL);CHKERRQ(ierr);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
230c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
231c4762a1bSJed Brown   appctx.gamma = .024;
232c4762a1bSJed Brown   appctx.kappa = .06;
233*60f0b76eSHong Zhang   appctx.aijpc = PETSC_FALSE;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* Create distributed array (DMDA) to manage parallel grid and vectors */
236c4762a1bSJed 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);
237c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
240c4762a1bSJed Brown   ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr);
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Extract global vectors from DMDA; then duplicate for remaining
243c4762a1bSJed Brown      vectors that are the same types */
244c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&appctx.U);CHKERRQ(ierr);
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* Create timestepping solver context */
247c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&appctx.ts);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = TSSetType(appctx.ts,TSCN);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = TSSetDM(appctx.ts,da);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = TSSetProblemType(appctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = TSSetEquationType(appctx.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
252c4762a1bSJed Brown   if (!implicitform) {
253c4762a1bSJed Brown     ierr = TSSetRHSFunction(appctx.ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr);
254c4762a1bSJed Brown     ierr = TSSetRHSJacobian(appctx.ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr);
255c4762a1bSJed Brown   } else {
256c4762a1bSJed Brown     ierr = TSSetIFunction(appctx.ts,NULL,IFunction,&appctx);CHKERRQ(ierr);
257c4762a1bSJed Brown     ierr = TSSetIJacobian(appctx.ts,NULL,NULL,IJacobian,&appctx);CHKERRQ(ierr);
258c4762a1bSJed Brown   }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* Set initial conditions */
261c4762a1bSJed Brown   ierr = InitialConditions(da,appctx.U);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = TSSetSolution(appctx.ts,appctx.U);CHKERRQ(ierr);
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* Set solver options */
265c4762a1bSJed Brown   ierr = TSSetMaxTime(appctx.ts,2000.0);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr = TSSetTimeStep(appctx.ts,0.5);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = TSSetExactFinalTime(appctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = TSSetFromOptions(appctx.ts);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   ierr = GenerateOBs(appctx.ts,appctx.U,&appctx);CHKERRQ(ierr);
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   if (!forwardonly) {
273c4762a1bSJed Brown     Tao           tao;
274c4762a1bSJed Brown     Vec           P;
275c4762a1bSJed Brown     Vec           lambda[1];
276956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
277c4762a1bSJed Brown     PetscLogStage opt_stage;
278956f8c0dSBarry Smith #endif
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     ierr = PetscLogStageRegister("Optimization",&opt_stage);CHKERRQ(ierr);
281c4762a1bSJed Brown     ierr = PetscLogStagePush(opt_stage);CHKERRQ(ierr);
282c4762a1bSJed Brown     if (perturbic == 1) {
283c4762a1bSJed Brown       ierr = PerturbedInitialConditions(da,appctx.U);CHKERRQ(ierr);
284c4762a1bSJed Brown     } else if (perturbic == 2) {
285c4762a1bSJed Brown       ierr = PerturbedInitialConditions2(da,appctx.U);CHKERRQ(ierr);
286c4762a1bSJed Brown     } else if (perturbic == 3) {
287c4762a1bSJed Brown       ierr = PerturbedInitialConditions3(da,appctx.U);CHKERRQ(ierr);
288c4762a1bSJed Brown     }
289c4762a1bSJed Brown 
290c4762a1bSJed Brown     ierr = VecDuplicate(appctx.U,&lambda[0]);CHKERRQ(ierr);
291c4762a1bSJed Brown     ierr = TSSetCostGradients(appctx.ts,1,lambda,NULL);CHKERRQ(ierr);
292c4762a1bSJed Brown 
293c4762a1bSJed Brown     /* Have the TS save its trajectory needed by TSAdjointSolve() */
294c4762a1bSJed Brown     ierr = TSSetSaveTrajectory(appctx.ts);CHKERRQ(ierr);
295c4762a1bSJed Brown 
296c4762a1bSJed Brown     /* Create TAO solver and set desired solution method */
297c4762a1bSJed Brown     ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
298c4762a1bSJed Brown     ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
299c4762a1bSJed Brown 
300c4762a1bSJed Brown     /* Set initial guess for TAO */
301c4762a1bSJed Brown     ierr = VecDuplicate(appctx.U,&P);CHKERRQ(ierr);
302c4762a1bSJed Brown     ierr = VecCopy(appctx.U,P);CHKERRQ(ierr);
303c4762a1bSJed Brown     ierr = TaoSetInitialVector(tao,P);CHKERRQ(ierr);
304c4762a1bSJed Brown 
305c4762a1bSJed Brown     /* Set routine for function and gradient evaluation */
306c4762a1bSJed Brown     ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionAndGradient,&appctx);CHKERRQ(ierr);
307c4762a1bSJed Brown 
308c4762a1bSJed Brown     /* Check for any TAO command line options */
309c4762a1bSJed Brown     ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
310c4762a1bSJed Brown 
311c4762a1bSJed Brown     ierr = TaoSolve(tao);CHKERRQ(ierr);
312c4762a1bSJed Brown     ierr = TaoDestroy(&tao);CHKERRQ(ierr);
313c4762a1bSJed Brown     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
314c4762a1bSJed Brown     ierr = VecDestroy(&P);CHKERRQ(ierr);
315c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
316c4762a1bSJed Brown   }
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /* Free work space.  All PETSc objects should be destroyed when they
319c4762a1bSJed Brown      are no longer needed. */
320c4762a1bSJed Brown   ierr = VecDestroy(&appctx.U);CHKERRQ(ierr);
321c4762a1bSJed Brown   ierr = TSDestroy(&appctx.ts);CHKERRQ(ierr);
322c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
323c4762a1bSJed Brown   ierr = PetscFinalize();
324c4762a1bSJed Brown   return ierr;
325c4762a1bSJed Brown }
326c4762a1bSJed Brown 
327c4762a1bSJed Brown /* ------------------------ TAO callbacks ---------------------------- */
328c4762a1bSJed Brown 
329c4762a1bSJed Brown /*
330c4762a1bSJed Brown    FormFunctionAndGradient - Evaluates the function and corresponding gradient.
331c4762a1bSJed Brown 
332c4762a1bSJed Brown    Input Parameters:
333c4762a1bSJed Brown    tao - the Tao context
334c4762a1bSJed Brown    P   - the input vector
335c4762a1bSJed Brown    ctx - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    Output Parameters:
338c4762a1bSJed Brown    f   - the newly evaluated function
339c4762a1bSJed Brown    G   - the newly evaluated gradient
340c4762a1bSJed Brown */
341c4762a1bSJed Brown PetscErrorCode FormFunctionAndGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
342c4762a1bSJed Brown {
343c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
344c4762a1bSJed Brown   PetscReal      soberr,timestep;
345c4762a1bSJed Brown   Vec            *lambda;
346c4762a1bSJed Brown   Vec            SDiff;
347c4762a1bSJed Brown   DM             da;
348c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN]="";
349c4762a1bSJed Brown   PetscViewer    viewer;
350c4762a1bSJed Brown   PetscErrorCode ierr;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   PetscFunctionBeginUser;
353c4762a1bSJed Brown   ierr = TSSetTime(appctx->ts,0.0);CHKERRQ(ierr);
354c4762a1bSJed Brown   ierr = TSGetTimeStep(appctx->ts,&timestep);CHKERRQ(ierr);
355c4762a1bSJed Brown   if (timestep<0) {
356c4762a1bSJed Brown     ierr = TSSetTimeStep(appctx->ts,-timestep);CHKERRQ(ierr);
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown   ierr = TSSetStepNumber(appctx->ts,0);CHKERRQ(ierr);
359c4762a1bSJed Brown   ierr = TSSetFromOptions(appctx->ts);CHKERRQ(ierr);
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   ierr = VecDuplicate(P,&SDiff);CHKERRQ(ierr);
362c4762a1bSJed Brown   ierr = VecCopy(P,appctx->U);CHKERRQ(ierr);
363c4762a1bSJed Brown   ierr = TSGetDM(appctx->ts,&da);CHKERRQ(ierr);
364c4762a1bSJed Brown   *f = 0;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   ierr = TSSolve(appctx->ts,appctx->U);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = PetscSNPrintf(filename,sizeof filename,"ex5opt.ob");CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = VecLoad(SDiff,viewer);CHKERRQ(ierr);
370c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
371c4762a1bSJed Brown   ierr = VecAYPX(SDiff,-1.,appctx->U);CHKERRQ(ierr);
372c4762a1bSJed Brown   ierr = VecDot(SDiff,SDiff,&soberr);CHKERRQ(ierr);
373c4762a1bSJed Brown   *f += soberr;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   ierr = TSGetCostGradients(appctx->ts,NULL,&lambda,NULL);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = VecSet(lambda[0],0.0);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = InitializeLambda(da,lambda[0],appctx->U,appctx);CHKERRQ(ierr);
378c4762a1bSJed Brown   ierr = TSAdjointSolve(appctx->ts);CHKERRQ(ierr);
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   ierr = VecCopy(lambda[0],G);CHKERRQ(ierr);
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   ierr = VecDestroy(&SDiff);CHKERRQ(ierr);
383c4762a1bSJed Brown   PetscFunctionReturn(0);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown /*TEST
387c4762a1bSJed Brown 
388c4762a1bSJed Brown    build:
389*60f0b76eSHong Zhang       depends: reaction_diffusion.c
390c4762a1bSJed Brown       requires: !complex !single
391c4762a1bSJed Brown 
392c4762a1bSJed Brown    test:
393c4762a1bSJed 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
394c4762a1bSJed Brown       output_file: output/ex5opt_ic_1.out
395c4762a1bSJed Brown 
396c4762a1bSJed Brown TEST*/
397