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