1c4762a1bSJed Brown static char help[] = "Demonstrates automatic, matrix-free Jacobian generation using ADOL-C for a time-dependent PDE in 2d, solved using implicit timestepping.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 5c4762a1bSJed Brown Concepts: Automatic differentiation using ADOL-C 6c4762a1bSJed Brown Concepts: Matrix-free methods 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 10c4762a1bSJed Brown 11c4762a1bSJed Brown For documentation on ADOL-C, see 12c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown /* ------------------------------------------------------------------------ 15c4762a1bSJed Brown See ../advection-diffusion-reaction/ex5 for a description of the problem 16c4762a1bSJed Brown ------------------------------------------------------------------------- */ 17c4762a1bSJed Brown 18c4762a1bSJed Brown #include <petscdmda.h> 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown #include "adolc-utils/init.cxx" 21c4762a1bSJed Brown #include "adolc-utils/matfree.cxx" 22c4762a1bSJed Brown #include <adolc/adolc.h> 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* (Passive) field for the two variables */ 25c4762a1bSJed Brown typedef struct { 26c4762a1bSJed Brown PetscScalar u,v; 27c4762a1bSJed Brown } Field; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Active field for the two variables */ 30c4762a1bSJed Brown typedef struct { 31c4762a1bSJed Brown adouble u,v; 32c4762a1bSJed Brown } AField; 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Application context */ 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown PetscReal D1,D2,gamma,kappa; 37c4762a1bSJed Brown AField **u_a,**f_a; 38c4762a1bSJed Brown AdolcCtx *adctx; /* Automatic differentation support */ 39c4762a1bSJed Brown } AppCtx; 40c4762a1bSJed Brown 41c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da,Vec U); 42c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y); 43c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr); 44c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr); 45c4762a1bSJed Brown extern PetscErrorCode IJacobianMatFree(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A_shell,Mat B,void *ctx); 46c4762a1bSJed Brown 47c4762a1bSJed Brown int main(int argc,char **argv) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown TS ts; /* ODE integrator */ 50c4762a1bSJed Brown Vec x,r; /* solution, residual */ 51c4762a1bSJed Brown DM da; 52c4762a1bSJed Brown AppCtx appctx; /* Application context */ 53c4762a1bSJed Brown AdolcMatCtx matctx; /* Matrix (free) context */ 54c4762a1bSJed Brown Vec lambda[1]; 55c4762a1bSJed Brown PetscBool forwardonly=PETSC_FALSE; 56c4762a1bSJed Brown Mat A; /* (Matrix free) Jacobian matrix */ 57c4762a1bSJed Brown PetscInt gxm,gym; 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 60c4762a1bSJed Brown Initialize program 61c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 62*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL)); 64c4762a1bSJed Brown PetscFunctionBeginUser; 65c4762a1bSJed Brown appctx.D1 = 8.0e-5; 66c4762a1bSJed Brown appctx.D2 = 4.0e-5; 67c4762a1bSJed Brown appctx.gamma = .024; 68c4762a1bSJed Brown appctx.kappa = .06; 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("df/dx forward",MAT_CLASSID,&matctx.event1)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("df/d(xdot) forward",MAT_CLASSID,&matctx.event2)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("df/dx reverse",MAT_CLASSID,&matctx.event3)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("df/d(xdot) reverse",MAT_CLASSID,&matctx.event4)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 76c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"u")); 815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,1,"v")); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 85c4762a1bSJed Brown vectors that are the same types 86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Create matrix free context and specify usage of PETSc-ADOL-C drivers 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATSHELL)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&A)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(A,&matctx)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))PetscAdolcIJacobianVectorProductIDMass)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void (*)(void))PetscAdolcIJacobianTransposeVectorProductIDMass)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&matctx.X)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&matctx.Xdot)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&matctx.localX0)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown Create timestepping solver context 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSCN)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)IFunctionLocalPassive,&appctx)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112c4762a1bSJed Brown Some data required for matrix-free context 113c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(da,NULL,NULL,NULL,&gxm,&gym,NULL)); 115c4762a1bSJed Brown matctx.m = 2*gxm*gym;matctx.n = 2*gxm*gym; /* Number of dependent and independent variables */ 116c4762a1bSJed Brown matctx.flg = PETSC_FALSE; /* Flag for reverse mode */ 117c4762a1bSJed Brown matctx.tag1 = 1; /* Tape identifier */ 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 120c4762a1bSJed Brown Trace function just once 121c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&appctx.adctx)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(IFunctionActive(ts,1.,x,matctx.Xdot,r,&appctx)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(appctx.adctx)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Set Jacobian. In this case, IJacobian simply acts to pass context 128c4762a1bSJed Brown information to the matrix-free Jacobian vector product. 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,A,A,IJacobianMatFree,&appctx)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Set initial conditions 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(da,x)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,x)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 140c4762a1bSJed Brown and set solver options 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142c4762a1bSJed Brown if (!forwardonly) { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(ts)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,200.0)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.5)); 146c4762a1bSJed Brown } else { 1475f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,2000.0)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,10)); 149c4762a1bSJed Brown } 1505f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 154c4762a1bSJed Brown Solve ODE system 155c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1565f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,x)); 157c4762a1bSJed Brown if (!forwardonly) { 158c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 159c4762a1bSJed Brown Start the Adjoint model 160c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&lambda[0])); 162c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 1635f80ce2aSJacob Faibussowitsch CHKERRQ(InitializeLambda(da,lambda[0],0.5,0.5)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(ts,1,lambda,NULL)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(ts)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&lambda[0])); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 171c4762a1bSJed Brown are no longer needed. 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&matctx.localX0)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&matctx.X)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&matctx.Xdot)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 181c4762a1bSJed Brown 182*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 183*b122ec5aSJacob Faibussowitsch return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U) 187c4762a1bSJed Brown { 188c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 189c4762a1bSJed Brown Field **u; 190c4762a1bSJed Brown PetscReal hx,hy,x,y; 191c4762a1bSJed Brown 192c4762a1bSJed Brown PetscFunctionBegin; 1935f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown hx = 2.5/(PetscReal)Mx; 196c4762a1bSJed Brown hy = 2.5/(PetscReal)My; 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* 199c4762a1bSJed Brown Get pointers to vector data 200c4762a1bSJed Brown */ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* 204c4762a1bSJed Brown Get local grid boundaries 205c4762a1bSJed Brown */ 2065f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown Compute function over the locally owned part of the grid 210c4762a1bSJed Brown */ 211c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 212c4762a1bSJed Brown y = j*hy; 213c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 214c4762a1bSJed Brown x = i*hx; 21566baab88SBarry Smith if (PetscApproximateGTE(x,1.0) && PetscApproximateLTE(x,1.5) && PetscApproximateGTE(y,1.0) && PetscApproximateLTE(y,1.5)) u[j][i].v = PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0)/4.0; 216c4762a1bSJed Brown else u[j][i].v = 0.0; 217c4762a1bSJed Brown 218c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0*u[j][i].v; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Restore vectors 224c4762a1bSJed Brown */ 2255f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 226c4762a1bSJed Brown PetscFunctionReturn(0); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y) 230c4762a1bSJed Brown { 231c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 232c4762a1bSJed Brown Field **l; 233c4762a1bSJed Brown 234410585f6SHong Zhang PetscFunctionBegin; 2355f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 236c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 237c4762a1bSJed Brown i = (PetscInt)(x*(Mx-1)); 238c4762a1bSJed Brown j = (PetscInt)(y*(My-1)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) { 242c4762a1bSJed Brown /* the i,j vertex is on this process */ 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,lambda,&l)); 244c4762a1bSJed Brown l[j][i].u = 1.0; 245c4762a1bSJed Brown l[j][i].v = 1.0; 2465f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,lambda,&l)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr) 252c4762a1bSJed Brown { 253c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 254c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym; 255c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 256c4762a1bSJed Brown PetscScalar uc,uxx,uyy,vc,vxx,vyy; 257c4762a1bSJed Brown 258c4762a1bSJed Brown PetscFunctionBegin; 259c4762a1bSJed Brown hx = 2.50/(PetscReal)(info->mx); sx = 1.0/(hx*hx); 260c4762a1bSJed Brown hy = 2.50/(PetscReal)(info->my); sy = 1.0/(hy*hy); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* Get local grid boundaries */ 263c4762a1bSJed Brown xs = info->xs; xm = info->xm; ys = info->ys; ym = info->ym; 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 266c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 267c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 268c4762a1bSJed Brown uc = u[j][i].u; 269c4762a1bSJed Brown uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 270c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 271c4762a1bSJed Brown vc = u[j][i].v; 272c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 273c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 274c4762a1bSJed Brown f[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc); 275c4762a1bSJed Brown f[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc; 276c4762a1bSJed Brown } 277c4762a1bSJed Brown } 2785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(16.0*xm*ym)); 279c4762a1bSJed Brown PetscFunctionReturn(0); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr) 283c4762a1bSJed Brown { 284c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 285c4762a1bSJed Brown DM da; 286c4762a1bSJed Brown DMDALocalInfo info; 287c4762a1bSJed Brown Field **u,**f,**udot; 288c4762a1bSJed Brown Vec localU; 289c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,gxs,gys,gxm,gym; 290c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 291c4762a1bSJed Brown adouble uc,uxx,uyy,vc,vxx,vyy; 292c4762a1bSJed Brown AField **f_a,*f_c,**u_a,*u_c; 293c4762a1bSJed Brown PetscScalar dummy; 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscFunctionBegin; 2965f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localU)); 299c4762a1bSJed Brown hx = 2.50/(PetscReal)(info.mx); sx = 1.0/(hx*hx); 300c4762a1bSJed Brown hy = 2.50/(PetscReal)(info.my); sy = 1.0/(hy*hy); 301c4762a1bSJed Brown xs = info.xs; xm = info.xm; gxs = info.gxs; gxm = info.gxm; 302c4762a1bSJed Brown ys = info.ys; ym = info.ym; gys = info.gys; gym = info.gym; 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* 305c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 306c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 307c4762a1bSJed Brown By placing code between these two statements, computations can be 308c4762a1bSJed Brown done while messages are in transition. 309c4762a1bSJed Brown */ 3105f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown Get pointers to vector data 315c4762a1bSJed Brown */ 3165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localU,&u)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,Udot,&udot)); 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Create contiguous 1-arrays of AFields 322c4762a1bSJed Brown 323c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 324c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 325c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 326c4762a1bSJed Brown */ 327c4762a1bSJed Brown u_c = new AField[info.gxm*info.gym]; 328c4762a1bSJed Brown f_c = new AField[info.gxm*info.gym]; 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 331c4762a1bSJed Brown u_a = new AField*[info.gym]; 332c4762a1bSJed Brown f_a = new AField*[info.gym]; 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 3355f80ce2aSJacob Faibussowitsch CHKERRQ(GiveGhostPoints(da,u_c,&u_a)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(GiveGhostPoints(da,f_c,&f_a)); 337c4762a1bSJed Brown 338c4762a1bSJed Brown trace_on(1); /* Start of active section on tape 1 */ 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* 341c4762a1bSJed Brown Mark independence 342c4762a1bSJed Brown 343c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 344c4762a1bSJed Brown other processors / on other boundaries. 345c4762a1bSJed Brown */ 346c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 347c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 348c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 349c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 350c4762a1bSJed Brown } 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 354c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 355c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 356c4762a1bSJed Brown uc = u_a[j][i].u; 357c4762a1bSJed Brown uxx = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx; 358c4762a1bSJed Brown uyy = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy; 359c4762a1bSJed Brown vc = u_a[j][i].v; 360c4762a1bSJed Brown vxx = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx; 361c4762a1bSJed Brown vyy = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy; 362c4762a1bSJed Brown f_a[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc); 363c4762a1bSJed Brown f_a[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc; 364c4762a1bSJed Brown } 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* 368c4762a1bSJed Brown Mark dependence 369c4762a1bSJed Brown 370c4762a1bSJed Brown NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming 371c4762a1bSJed Brown the Jacobian later. 372c4762a1bSJed Brown */ 373c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 374c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 375c4762a1bSJed Brown if ((i < xs) || (i >= xs+xm) || (j < ys) || (j >= ys+ym)) { 376c4762a1bSJed Brown f_a[j][i].u >>= dummy; 377c4762a1bSJed Brown f_a[j][i].v >>= dummy; 378c4762a1bSJed Brown } else { 379c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 380c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 381c4762a1bSJed Brown } 382c4762a1bSJed Brown } 383c4762a1bSJed Brown } 384c4762a1bSJed Brown trace_off(); /* End of active section */ 3855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(16.0*xm*ym)); 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* Restore vectors */ 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,Udot,&udot)); 391c4762a1bSJed Brown 3925f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localU)); 393410585f6SHong Zhang 394c4762a1bSJed Brown /* Destroy AFields appropriately */ 395c4762a1bSJed Brown f_a += info.gys; 396c4762a1bSJed Brown u_a += info.gys; 397c4762a1bSJed Brown delete[] f_a; 398c4762a1bSJed Brown delete[] u_a; 399c4762a1bSJed Brown delete[] f_c; 400c4762a1bSJed Brown delete[] u_c; 401c4762a1bSJed Brown PetscFunctionReturn(0); 402c4762a1bSJed Brown } 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* 405c4762a1bSJed Brown Simply acts to pass TS information to the AdolcMatCtx 406c4762a1bSJed Brown */ 407c4762a1bSJed Brown PetscErrorCode IJacobianMatFree(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A_shell,Mat B,void *ctx) 408c4762a1bSJed Brown { 409c4762a1bSJed Brown AdolcMatCtx *mctx; 410c4762a1bSJed Brown DM da; 411c4762a1bSJed Brown 412c4762a1bSJed Brown PetscFunctionBeginUser; 4135f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A_shell,&mctx)); 414c4762a1bSJed Brown 415c4762a1bSJed Brown mctx->time = t; 416c4762a1bSJed Brown mctx->shift = a; 417c4762a1bSJed Brown if (mctx->ts != ts) mctx->ts = ts; 4185f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(X,mctx->X)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(Xdot,mctx->Xdot)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,mctx->X,INSERT_VALUES,mctx->localX0)); 4225f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,mctx->X,INSERT_VALUES,mctx->localX0)); 423c4762a1bSJed Brown PetscFunctionReturn(0); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown 426c4762a1bSJed Brown /*TEST 427c4762a1bSJed Brown 428c4762a1bSJed Brown build: 429c4762a1bSJed Brown requires: double !complex adolc 430c4762a1bSJed Brown 431c4762a1bSJed Brown test: 432c4762a1bSJed Brown suffix: 1 433c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian 434c4762a1bSJed Brown output_file: output/adr_ex5adj_mf_1.out 435c4762a1bSJed Brown 436c4762a1bSJed Brown test: 437c4762a1bSJed Brown suffix: 2 438c4762a1bSJed Brown nsize: 4 439c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor 440c4762a1bSJed Brown output_file: output/adr_ex5adj_mf_2.out 441c4762a1bSJed Brown 442c4762a1bSJed Brown TEST*/ 443