1*c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 5*c4762a1bSJed Brown Concepts: Automatic differentiation using ADOL-C 6*c4762a1bSJed Brown */ 7*c4762a1bSJed Brown /* 8*c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown For documentation on ADOL-C, see 11*c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-colpack 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown For documentation on ColPack, see 16*c4762a1bSJed Brown $PETSC_ARCH/externalpackages/git.colpack/README.md 17*c4762a1bSJed Brown */ 18*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 19*c4762a1bSJed Brown See ../advection-diffusion-reaction/ex5 for a description of the problem 20*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown /* 23*c4762a1bSJed Brown Runtime options: 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown Solver: 26*c4762a1bSJed Brown -forwardonly - Run the forward simulation without adjoint. 27*c4762a1bSJed Brown -implicitform - Provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used. 28*c4762a1bSJed Brown -aijpc - Set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL). 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown Jacobian generation: 31*c4762a1bSJed Brown -jacobian_by_hand - Use the hand-coded Jacobian of ex5.c, rather than generating it automatically. 32*c4762a1bSJed Brown -no_annotation - Do not annotate ADOL-C active variables. (Should be used alongside -jacobian_by_hand.) 33*c4762a1bSJed Brown -adolc_sparse - Calculate Jacobian in compressed format, using ADOL-C sparse functionality. 34*c4762a1bSJed Brown -adolc_sparse_view - Print sparsity pattern used by -adolc_sparse option. 35*c4762a1bSJed Brown */ 36*c4762a1bSJed Brown /* 37*c4762a1bSJed Brown NOTE: If -adolc_sparse option is used, at least four processors should be used, so that no processor owns boundaries which are 38*c4762a1bSJed Brown identified by the periodic boundary conditions. The number of grid points in both x- and y-directions should be multiples 39*c4762a1bSJed Brown of 5, in order for the 5-point stencil to be cleanly parallelised. 40*c4762a1bSJed Brown */ 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown #include <petscdmda.h> 43*c4762a1bSJed Brown #include <petscts.h> 44*c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 45*c4762a1bSJed Brown #include <adolc/adolc.h> 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown /* (Passive) field for the two variables */ 48*c4762a1bSJed Brown typedef struct { 49*c4762a1bSJed Brown PetscScalar u,v; 50*c4762a1bSJed Brown } Field; 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown /* Active field for the two variables */ 53*c4762a1bSJed Brown typedef struct { 54*c4762a1bSJed Brown adouble u,v; 55*c4762a1bSJed Brown } AField; 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* Application context */ 58*c4762a1bSJed Brown typedef struct { 59*c4762a1bSJed Brown PetscReal D1,D2,gamma,kappa; 60*c4762a1bSJed Brown AField **u_a,**f_a; 61*c4762a1bSJed Brown PetscBool aijpc; 62*c4762a1bSJed Brown AdolcCtx *adctx; /* Automatic differentation support */ 63*c4762a1bSJed Brown } AppCtx; 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da,Vec U); 66*c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y); 67*c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr); 68*c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr); 69*c4762a1bSJed Brown extern PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr); 70*c4762a1bSJed Brown extern PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr); 71*c4762a1bSJed Brown extern PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx); 72*c4762a1bSJed Brown extern PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx); 73*c4762a1bSJed Brown extern PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx); 74*c4762a1bSJed Brown extern PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown int main(int argc,char **argv) 77*c4762a1bSJed Brown { 78*c4762a1bSJed Brown TS ts; /* ODE integrator */ 79*c4762a1bSJed Brown Vec x,r,xdot; /* solution, residual, derivative */ 80*c4762a1bSJed Brown PetscErrorCode ierr; 81*c4762a1bSJed Brown DM da; 82*c4762a1bSJed Brown AppCtx appctx; 83*c4762a1bSJed Brown AdolcCtx *adctx; 84*c4762a1bSJed Brown Vec lambda[1]; 85*c4762a1bSJed Brown PetscBool forwardonly=PETSC_FALSE,implicitform=PETSC_FALSE,byhand=PETSC_FALSE; 86*c4762a1bSJed Brown PetscInt gxm,gym,i,dofs = 2,ctrl[3] = {0,0,0}; 87*c4762a1bSJed Brown PetscScalar **Seed = NULL,**Rec = NULL,*u_vec; 88*c4762a1bSJed Brown unsigned int **JP = NULL; 89*c4762a1bSJed Brown ISColoring iscoloring; 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 92*c4762a1bSJed Brown ierr = PetscNew(&adctx);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr); 95*c4762a1bSJed Brown appctx.aijpc = PETSC_FALSE,adctx->no_an = PETSC_FALSE,adctx->sparse = PETSC_FALSE,adctx->sparse_view = PETSC_FALSE;adctx->sparse_view_done = PETSC_FALSE; 96*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-no_annotation",&adctx->no_an,NULL);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-jacobian_by_hand",&byhand,NULL);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-adolc_sparse",&adctx->sparse,NULL);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-adolc_sparse_view",&adctx->sparse_view,NULL);CHKERRQ(ierr); 101*c4762a1bSJed Brown appctx.D1 = 8.0e-5; 102*c4762a1bSJed Brown appctx.D2 = 4.0e-5; 103*c4762a1bSJed Brown appctx.gamma = .024; 104*c4762a1bSJed Brown appctx.kappa = .06; 105*c4762a1bSJed Brown appctx.adctx = adctx; 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108*c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 109*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr); 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117*c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 118*c4762a1bSJed Brown vectors that are the same types 119*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecDuplicate(x,&xdot);CHKERRQ(ierr); 123*c4762a1bSJed Brown 124*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125*c4762a1bSJed Brown Create timestepping solver context 126*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 131*c4762a1bSJed Brown if (!implicitform) { 132*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&appctx);CHKERRQ(ierr); 133*c4762a1bSJed Brown } else { 134*c4762a1bSJed Brown ierr = DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)IFunctionLocalPassive,&appctx);CHKERRQ(ierr); 135*c4762a1bSJed Brown } 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown if (!adctx->no_an) { 138*c4762a1bSJed Brown ierr = DMDAGetGhostCorners(da,NULL,NULL,NULL,&gxm,&gym,NULL);CHKERRQ(ierr); 139*c4762a1bSJed Brown adctx->m = dofs*gxm*gym;adctx->n = dofs*gxm*gym; /* Number of dependent and independent variables */ 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142*c4762a1bSJed Brown Trace function(s) just once 143*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144*c4762a1bSJed Brown ierr = PetscMalloc1(adctx->n,&u_vec);CHKERRQ(ierr); 145*c4762a1bSJed Brown if (!implicitform) { 146*c4762a1bSJed Brown ierr = RHSFunctionActive(ts,1.0,x,r,&appctx);CHKERRQ(ierr); 147*c4762a1bSJed Brown } else { 148*c4762a1bSJed Brown ierr = IFunctionActive(ts,1.0,x,xdot,r,&appctx);CHKERRQ(ierr); 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown 151*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152*c4762a1bSJed Brown In the case where ADOL-C generates the Jacobian in compressed format, 153*c4762a1bSJed Brown seed and recovery matrices are required. Since the sparsity structure 154*c4762a1bSJed Brown of the Jacobian does not change over the course of the time 155*c4762a1bSJed Brown integration, we can save computational effort by only generating 156*c4762a1bSJed Brown these objects once. 157*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158*c4762a1bSJed Brown if (adctx->sparse) { 159*c4762a1bSJed Brown 160*c4762a1bSJed Brown /* 161*c4762a1bSJed Brown Generate sparsity pattern 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown One way to do this is to use the Jacobian pattern driver `jac_pat` 164*c4762a1bSJed Brown provided by ColPack. Otherwise, it is the responsibility of the user 165*c4762a1bSJed Brown to obtain the sparsity pattern. 166*c4762a1bSJed Brown */ 167*c4762a1bSJed Brown ierr = PetscMalloc1(adctx->n,&u_vec);CHKERRQ(ierr); 168*c4762a1bSJed Brown JP = (unsigned int **) malloc(adctx->m*sizeof(unsigned int*)); 169*c4762a1bSJed Brown jac_pat(1,adctx->m,adctx->n,u_vec,JP,ctrl); 170*c4762a1bSJed Brown if (adctx->sparse_view) { 171*c4762a1bSJed Brown ierr = PrintSparsity(MPI_COMM_WORLD,adctx->m,JP);CHKERRQ(ierr); 172*c4762a1bSJed Brown } 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown /* 175*c4762a1bSJed Brown Extract a column colouring 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown For problems using DMDA, colourings can be extracted directly, as 178*c4762a1bSJed Brown follows. There are also ColPack drivers available. Otherwise, it is the 179*c4762a1bSJed Brown responsibility of the user to obtain a suitable colouring. 180*c4762a1bSJed Brown */ 181*c4762a1bSJed Brown ierr = DMCreateColoring(da,IS_COLORING_LOCAL,&iscoloring);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&adctx->p,NULL);CHKERRQ(ierr); 183*c4762a1bSJed Brown 184*c4762a1bSJed Brown /* Generate seed matrix to propagate through the forward mode of AD */ 185*c4762a1bSJed Brown ierr = AdolcMalloc2(adctx->n,adctx->p,&Seed);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = GenerateSeedMatrix(iscoloring,Seed);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* 190*c4762a1bSJed Brown Generate recovery matrix, which is used to recover the Jacobian from 191*c4762a1bSJed Brown compressed format */ 192*c4762a1bSJed Brown ierr = AdolcMalloc2(adctx->m,adctx->p,&Rec);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = GetRecoveryMatrix(Seed,JP,adctx->m,adctx->p,Rec);CHKERRQ(ierr); 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown /* Store results and free workspace */ 196*c4762a1bSJed Brown adctx->Rec = Rec; 197*c4762a1bSJed Brown for (i=0;i<adctx->m;i++) 198*c4762a1bSJed Brown free(JP[i]); 199*c4762a1bSJed Brown free(JP); 200*c4762a1bSJed Brown ierr = PetscFree(u_vec);CHKERRQ(ierr); 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown } else { 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown /* 205*c4762a1bSJed Brown In 'full' Jacobian mode, propagate an identity matrix through the 206*c4762a1bSJed Brown forward mode of AD. 207*c4762a1bSJed Brown */ 208*c4762a1bSJed Brown adctx->p = adctx->n; 209*c4762a1bSJed Brown ierr = AdolcMalloc2(adctx->n,adctx->p,&Seed);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = Identity(adctx->n,Seed);CHKERRQ(ierr); 211*c4762a1bSJed Brown } 212*c4762a1bSJed Brown adctx->Seed = Seed; 213*c4762a1bSJed Brown } 214*c4762a1bSJed Brown 215*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 216*c4762a1bSJed Brown Set Jacobian 217*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 218*c4762a1bSJed Brown if (!implicitform) { 219*c4762a1bSJed Brown if (!byhand) { 220*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianAdolc,&appctx);CHKERRQ(ierr); 221*c4762a1bSJed Brown } else { 222*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianByHand,&appctx);CHKERRQ(ierr); 223*c4762a1bSJed Brown } 224*c4762a1bSJed Brown } else { 225*c4762a1bSJed Brown if (appctx.aijpc) { 226*c4762a1bSJed Brown Mat A,B; 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown ierr = DMSetMatType(da,MATSELL);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 231*c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */ 232*c4762a1bSJed Brown if (!byhand) { 233*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,B,IJacobianAdolc,&appctx);CHKERRQ(ierr); 234*c4762a1bSJed Brown } else { 235*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,B,IJacobianByHand,&appctx);CHKERRQ(ierr); 236*c4762a1bSJed Brown } 237*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 239*c4762a1bSJed Brown } else { 240*c4762a1bSJed Brown if (!byhand) { 241*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,NULL,NULL,IJacobianAdolc,&appctx);CHKERRQ(ierr); 242*c4762a1bSJed Brown } else { 243*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,NULL,NULL,IJacobianByHand,&appctx);CHKERRQ(ierr); 244*c4762a1bSJed Brown } 245*c4762a1bSJed Brown } 246*c4762a1bSJed Brown } 247*c4762a1bSJed Brown 248*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249*c4762a1bSJed Brown Set initial conditions 250*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 251*c4762a1bSJed Brown ierr = InitialConditions(da,x);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 253*c4762a1bSJed Brown 254*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 255*c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 256*c4762a1bSJed Brown and set solver options 257*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 258*c4762a1bSJed Brown if (!forwardonly) { 259*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 260*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,200.0);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.5);CHKERRQ(ierr); 262*c4762a1bSJed Brown } else { 263*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,2000.0);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,10);CHKERRQ(ierr); 265*c4762a1bSJed Brown } 266*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 270*c4762a1bSJed Brown Solve ODE system 271*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 272*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 273*c4762a1bSJed Brown if (!forwardonly) { 274*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 275*c4762a1bSJed Brown Start the Adjoint model 276*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 277*c4762a1bSJed Brown ierr = VecDuplicate(x,&lambda[0]);CHKERRQ(ierr); 278*c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 279*c4762a1bSJed Brown ierr = InitializeLambda(da,lambda[0],0.5,0.5);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,lambda,NULL);CHKERRQ(ierr); 281*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 282*c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 283*c4762a1bSJed Brown } 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 286*c4762a1bSJed Brown Free work space. 287*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 288*c4762a1bSJed Brown ierr = VecDestroy(&xdot);CHKERRQ(ierr); 289*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 290*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 291*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 292*c4762a1bSJed Brown if (!adctx->no_an) { 293*c4762a1bSJed Brown if (adctx->sparse) 294*c4762a1bSJed Brown ierr = AdolcFree2(Rec);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = AdolcFree2(Seed);CHKERRQ(ierr); 296*c4762a1bSJed Brown } 297*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = PetscFree(adctx);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = PetscFinalize(); 300*c4762a1bSJed Brown return ierr; 301*c4762a1bSJed Brown } 302*c4762a1bSJed Brown 303*c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U) 304*c4762a1bSJed Brown { 305*c4762a1bSJed Brown PetscErrorCode ierr; 306*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 307*c4762a1bSJed Brown Field **u; 308*c4762a1bSJed Brown PetscReal hx,hy,x,y; 309*c4762a1bSJed Brown 310*c4762a1bSJed Brown PetscFunctionBegin; 311*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); 312*c4762a1bSJed Brown 313*c4762a1bSJed Brown hx = 2.5/(PetscReal)Mx; 314*c4762a1bSJed Brown hy = 2.5/(PetscReal)My; 315*c4762a1bSJed Brown 316*c4762a1bSJed Brown /* 317*c4762a1bSJed Brown Get pointers to vector data 318*c4762a1bSJed Brown */ 319*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown /* 322*c4762a1bSJed Brown Get local grid boundaries 323*c4762a1bSJed Brown */ 324*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 325*c4762a1bSJed Brown 326*c4762a1bSJed Brown /* 327*c4762a1bSJed Brown Compute function over the locally owned part of the grid 328*c4762a1bSJed Brown */ 329*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 330*c4762a1bSJed Brown y = j*hy; 331*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 332*c4762a1bSJed Brown x = i*hx; 333*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); 334*c4762a1bSJed Brown else u[j][i].v = 0.0; 335*c4762a1bSJed Brown 336*c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0*u[j][i].v; 337*c4762a1bSJed Brown } 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown /* 341*c4762a1bSJed Brown Restore vectors 342*c4762a1bSJed Brown */ 343*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 344*c4762a1bSJed Brown PetscFunctionReturn(0); 345*c4762a1bSJed Brown } 346*c4762a1bSJed Brown 347*c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y) 348*c4762a1bSJed Brown { 349*c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 350*c4762a1bSJed Brown PetscErrorCode ierr; 351*c4762a1bSJed Brown Field **l; 352*c4762a1bSJed Brown PetscFunctionBegin; 353*c4762a1bSJed Brown 354*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); 355*c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 356*c4762a1bSJed Brown i = (PetscInt)(x*(Mx-1)); 357*c4762a1bSJed Brown j = (PetscInt)(y*(My-1)); 358*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 359*c4762a1bSJed Brown 360*c4762a1bSJed Brown if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) { 361*c4762a1bSJed Brown /* the i,j vertex is on this process */ 362*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,lambda,&l);CHKERRQ(ierr); 363*c4762a1bSJed Brown l[j][i].u = 1.0; 364*c4762a1bSJed Brown l[j][i].v = 1.0; 365*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,lambda,&l);CHKERRQ(ierr); 366*c4762a1bSJed Brown } 367*c4762a1bSJed Brown PetscFunctionReturn(0); 368*c4762a1bSJed Brown } 369*c4762a1bSJed Brown 370*c4762a1bSJed Brown PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr) 371*c4762a1bSJed Brown { 372*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 373*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym; 374*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 375*c4762a1bSJed Brown PetscScalar uc,uxx,uyy,vc,vxx,vyy; 376*c4762a1bSJed Brown PetscErrorCode ierr; 377*c4762a1bSJed Brown 378*c4762a1bSJed Brown PetscFunctionBegin; 379*c4762a1bSJed Brown hx = 2.50/(PetscReal)(info->mx); sx = 1.0/(hx*hx); 380*c4762a1bSJed Brown hy = 2.50/(PetscReal)(info->my); sy = 1.0/(hy*hy); 381*c4762a1bSJed Brown 382*c4762a1bSJed Brown /* Get local grid boundaries */ 383*c4762a1bSJed Brown xs = info->xs; xm = info->xm; ys = info->ys; ym = info->ym; 384*c4762a1bSJed Brown 385*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 386*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 387*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 388*c4762a1bSJed Brown uc = u[j][i].u; 389*c4762a1bSJed Brown uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 390*c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 391*c4762a1bSJed Brown vc = u[j][i].v; 392*c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 393*c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 394*c4762a1bSJed Brown f[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc); 395*c4762a1bSJed Brown f[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc; 396*c4762a1bSJed Brown } 397*c4762a1bSJed Brown } 398*c4762a1bSJed Brown ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr); 399*c4762a1bSJed Brown PetscFunctionReturn(0); 400*c4762a1bSJed Brown } 401*c4762a1bSJed Brown 402*c4762a1bSJed Brown PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr) 403*c4762a1bSJed Brown { 404*c4762a1bSJed Brown PetscErrorCode ierr; 405*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 406*c4762a1bSJed Brown DM da; 407*c4762a1bSJed Brown DMDALocalInfo info; 408*c4762a1bSJed Brown Field **u,**f,**udot; 409*c4762a1bSJed Brown Vec localU; 410*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,gxs,gys,gxm,gym; 411*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 412*c4762a1bSJed Brown adouble uc,uxx,uyy,vc,vxx,vyy; 413*c4762a1bSJed Brown AField **f_a,*f_c,**u_a,*u_c; 414*c4762a1bSJed Brown PetscScalar dummy; 415*c4762a1bSJed Brown 416*c4762a1bSJed Brown PetscFunctionBegin; 417*c4762a1bSJed Brown 418*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 419*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 420*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 421*c4762a1bSJed Brown hx = 2.50/(PetscReal)(info.mx); sx = 1.0/(hx*hx); 422*c4762a1bSJed Brown hy = 2.50/(PetscReal)(info.my); sy = 1.0/(hy*hy); 423*c4762a1bSJed Brown xs = info.xs; xm = info.xm; gxs = info.gxs; gxm = info.gxm; 424*c4762a1bSJed Brown ys = info.ys; ym = info.ym; gys = info.gys; gym = info.gym; 425*c4762a1bSJed Brown 426*c4762a1bSJed Brown /* 427*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 428*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 429*c4762a1bSJed Brown By placing code between these two statements, computations can be 430*c4762a1bSJed Brown done while messages are in transition. 431*c4762a1bSJed Brown */ 432*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 433*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 434*c4762a1bSJed Brown 435*c4762a1bSJed Brown /* 436*c4762a1bSJed Brown Get pointers to vector data 437*c4762a1bSJed Brown */ 438*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 439*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 440*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr); 441*c4762a1bSJed Brown 442*c4762a1bSJed Brown /* 443*c4762a1bSJed Brown Create contiguous 1-arrays of AFields 444*c4762a1bSJed Brown 445*c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 446*c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 447*c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 448*c4762a1bSJed Brown */ 449*c4762a1bSJed Brown u_c = new AField[info.gxm*info.gym]; 450*c4762a1bSJed Brown f_c = new AField[info.gxm*info.gym]; 451*c4762a1bSJed Brown 452*c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 453*c4762a1bSJed Brown u_a = new AField*[info.gym]; 454*c4762a1bSJed Brown f_a = new AField*[info.gym]; 455*c4762a1bSJed Brown 456*c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 457*c4762a1bSJed Brown ierr = GiveGhostPoints(da,u_c,&u_a);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = GiveGhostPoints(da,f_c,&f_a);CHKERRQ(ierr); 459*c4762a1bSJed Brown 460*c4762a1bSJed Brown trace_on(1); /* Start of active section on tape 1 */ 461*c4762a1bSJed Brown 462*c4762a1bSJed Brown /* 463*c4762a1bSJed Brown Mark independence 464*c4762a1bSJed Brown 465*c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 466*c4762a1bSJed Brown other processors / on other boundaries. 467*c4762a1bSJed Brown */ 468*c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 469*c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 470*c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 471*c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 472*c4762a1bSJed Brown } 473*c4762a1bSJed Brown } 474*c4762a1bSJed Brown 475*c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 476*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 477*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 478*c4762a1bSJed Brown uc = u_a[j][i].u; 479*c4762a1bSJed Brown uxx = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx; 480*c4762a1bSJed Brown uyy = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy; 481*c4762a1bSJed Brown vc = u_a[j][i].v; 482*c4762a1bSJed Brown vxx = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx; 483*c4762a1bSJed Brown vyy = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy; 484*c4762a1bSJed Brown f_a[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc); 485*c4762a1bSJed Brown f_a[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc; 486*c4762a1bSJed Brown } 487*c4762a1bSJed Brown } 488*c4762a1bSJed Brown 489*c4762a1bSJed Brown /* 490*c4762a1bSJed Brown Mark dependence 491*c4762a1bSJed Brown 492*c4762a1bSJed Brown NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming 493*c4762a1bSJed Brown the Jacobian later. 494*c4762a1bSJed Brown */ 495*c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 496*c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 497*c4762a1bSJed Brown if ((i < xs) || (i >= xs+xm) || (j < ys) || (j >= ys+ym)) { 498*c4762a1bSJed Brown f_a[j][i].u >>= dummy; 499*c4762a1bSJed Brown f_a[j][i].v >>= dummy; 500*c4762a1bSJed Brown } else { 501*c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 502*c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 503*c4762a1bSJed Brown } 504*c4762a1bSJed Brown } 505*c4762a1bSJed Brown } 506*c4762a1bSJed Brown trace_off(); /* End of active section */ 507*c4762a1bSJed Brown ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr); 508*c4762a1bSJed Brown 509*c4762a1bSJed Brown /* Restore vectors */ 510*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 511*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 512*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr); 513*c4762a1bSJed Brown 514*c4762a1bSJed Brown /* Destroy AFields appropriately */ 515*c4762a1bSJed Brown f_a += info.gys; 516*c4762a1bSJed Brown u_a += info.gys; 517*c4762a1bSJed Brown delete[] f_a; 518*c4762a1bSJed Brown delete[] u_a; 519*c4762a1bSJed Brown delete[] f_c; 520*c4762a1bSJed Brown delete[] u_c; 521*c4762a1bSJed Brown 522*c4762a1bSJed Brown PetscFunctionReturn(0); 523*c4762a1bSJed Brown } 524*c4762a1bSJed Brown 525*c4762a1bSJed Brown PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr) 526*c4762a1bSJed Brown { 527*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 528*c4762a1bSJed Brown DM da; 529*c4762a1bSJed Brown PetscErrorCode ierr; 530*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 531*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 532*c4762a1bSJed Brown PetscScalar uc,uxx,uyy,vc,vxx,vyy; 533*c4762a1bSJed Brown Field **u,**f; 534*c4762a1bSJed Brown Vec localU,localF; 535*c4762a1bSJed Brown 536*c4762a1bSJed Brown PetscFunctionBegin; 537*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 538*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); 539*c4762a1bSJed Brown hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx); 540*c4762a1bSJed Brown hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy); 541*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 542*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr); 543*c4762a1bSJed Brown 544*c4762a1bSJed Brown /* 545*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 546*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 547*c4762a1bSJed Brown By placing code between these two statements, computations can be 548*c4762a1bSJed Brown done while messages are in transition. 549*c4762a1bSJed Brown */ 550*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 551*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 552*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); // NOTE (1): See (2) below 553*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF);CHKERRQ(ierr); 554*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF);CHKERRQ(ierr); 555*c4762a1bSJed Brown 556*c4762a1bSJed Brown /* 557*c4762a1bSJed Brown Get pointers to vector data 558*c4762a1bSJed Brown */ 559*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 560*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,localF,&f);CHKERRQ(ierr); 561*c4762a1bSJed Brown 562*c4762a1bSJed Brown /* 563*c4762a1bSJed Brown Get local grid boundaries 564*c4762a1bSJed Brown */ 565*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 566*c4762a1bSJed Brown 567*c4762a1bSJed Brown /* 568*c4762a1bSJed Brown Compute function over the locally owned part of the grid 569*c4762a1bSJed Brown */ 570*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 571*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 572*c4762a1bSJed Brown uc = u[j][i].u; 573*c4762a1bSJed Brown uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 574*c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 575*c4762a1bSJed Brown vc = u[j][i].v; 576*c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 577*c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 578*c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc); 579*c4762a1bSJed Brown f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; 580*c4762a1bSJed Brown } 581*c4762a1bSJed Brown } 582*c4762a1bSJed Brown 583*c4762a1bSJed Brown /* 584*c4762a1bSJed Brown Gather global vector, using the 2-step process 585*c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd(). 586*c4762a1bSJed Brown 587*c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or 588*c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above). 589*c4762a1bSJed Brown */ 590*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr); 591*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr); 592*c4762a1bSJed Brown 593*c4762a1bSJed Brown /* 594*c4762a1bSJed Brown Restore vectors 595*c4762a1bSJed Brown */ 596*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,localF,&f);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 598*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr); 599*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr); 601*c4762a1bSJed Brown PetscFunctionReturn(0); 602*c4762a1bSJed Brown } 603*c4762a1bSJed Brown 604*c4762a1bSJed Brown PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr) 605*c4762a1bSJed Brown { 606*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 607*c4762a1bSJed Brown DM da; 608*c4762a1bSJed Brown PetscErrorCode ierr; 609*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,gxs,gys,gxm,gym,Mx,My; 610*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 611*c4762a1bSJed Brown AField **f_a,*f_c,**u_a,*u_c; 612*c4762a1bSJed Brown adouble uc,uxx,uyy,vc,vxx,vyy; 613*c4762a1bSJed Brown Field **u,**f; 614*c4762a1bSJed Brown Vec localU,localF; 615*c4762a1bSJed Brown 616*c4762a1bSJed Brown PetscFunctionBegin; 617*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 618*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); 619*c4762a1bSJed Brown hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx); 620*c4762a1bSJed Brown hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy); 621*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 622*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr); 623*c4762a1bSJed Brown 624*c4762a1bSJed Brown /* 625*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 626*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 627*c4762a1bSJed Brown By placing code between these two statements, computations can be 628*c4762a1bSJed Brown done while messages are in transition. 629*c4762a1bSJed Brown */ 630*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 631*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 632*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); // NOTE (1): See (2) below 633*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF);CHKERRQ(ierr); 634*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF);CHKERRQ(ierr); 635*c4762a1bSJed Brown 636*c4762a1bSJed Brown /* 637*c4762a1bSJed Brown Get pointers to vector data 638*c4762a1bSJed Brown */ 639*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 640*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,localF,&f);CHKERRQ(ierr); 641*c4762a1bSJed Brown 642*c4762a1bSJed Brown /* 643*c4762a1bSJed Brown Get local and ghosted grid boundaries 644*c4762a1bSJed Brown */ 645*c4762a1bSJed Brown ierr = DMDAGetGhostCorners(da,&gxs,&gys,NULL,&gxm,&gym,NULL);CHKERRQ(ierr); 646*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 647*c4762a1bSJed Brown 648*c4762a1bSJed Brown /* 649*c4762a1bSJed Brown Create contiguous 1-arrays of AFields 650*c4762a1bSJed Brown 651*c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 652*c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 653*c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 654*c4762a1bSJed Brown */ 655*c4762a1bSJed Brown u_c = new AField[gxm*gym]; 656*c4762a1bSJed Brown f_c = new AField[gxm*gym]; 657*c4762a1bSJed Brown 658*c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 659*c4762a1bSJed Brown u_a = new AField*[gym]; 660*c4762a1bSJed Brown f_a = new AField*[gym]; 661*c4762a1bSJed Brown 662*c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 663*c4762a1bSJed Brown ierr = GiveGhostPoints(da,u_c,&u_a);CHKERRQ(ierr); 664*c4762a1bSJed Brown ierr = GiveGhostPoints(da,f_c,&f_a);CHKERRQ(ierr); 665*c4762a1bSJed Brown 666*c4762a1bSJed Brown /* 667*c4762a1bSJed Brown Compute function over the locally owned part of the grid 668*c4762a1bSJed Brown */ 669*c4762a1bSJed Brown trace_on(1); // ----------------------------------------------- Start of active section 670*c4762a1bSJed Brown 671*c4762a1bSJed Brown /* 672*c4762a1bSJed Brown Mark independence 673*c4762a1bSJed Brown 674*c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 675*c4762a1bSJed Brown other processors / on other boundaries. 676*c4762a1bSJed Brown */ 677*c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 678*c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 679*c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 680*c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 681*c4762a1bSJed Brown } 682*c4762a1bSJed Brown } 683*c4762a1bSJed Brown 684*c4762a1bSJed Brown /* 685*c4762a1bSJed Brown Compute function over the locally owned part of the grid 686*c4762a1bSJed Brown */ 687*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 688*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 689*c4762a1bSJed Brown uc = u_a[j][i].u; 690*c4762a1bSJed Brown uxx = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx; 691*c4762a1bSJed Brown uyy = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy; 692*c4762a1bSJed Brown vc = u_a[j][i].v; 693*c4762a1bSJed Brown vxx = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx; 694*c4762a1bSJed Brown vyy = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy; 695*c4762a1bSJed Brown f_a[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc); 696*c4762a1bSJed Brown f_a[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; 697*c4762a1bSJed Brown } 698*c4762a1bSJed Brown } 699*c4762a1bSJed Brown /* 700*c4762a1bSJed Brown Mark dependence 701*c4762a1bSJed Brown 702*c4762a1bSJed Brown NOTE: Ghost points are marked as dependent in order to vastly simplify index notation 703*c4762a1bSJed Brown during Jacobian assembly. 704*c4762a1bSJed Brown */ 705*c4762a1bSJed Brown for (j=gys; j<gys+gym; j++) { 706*c4762a1bSJed Brown for (i=gxs; i<gxs+gxm; i++) { 707*c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 708*c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 709*c4762a1bSJed Brown } 710*c4762a1bSJed Brown } 711*c4762a1bSJed Brown trace_off(); // ----------------------------------------------- End of active section 712*c4762a1bSJed Brown 713*c4762a1bSJed Brown /* Test zeroth order scalar evaluation in ADOL-C gives the same result */ 714*c4762a1bSJed Brown // if (appctx->adctx->zos) { 715*c4762a1bSJed Brown // ierr = TestZOS2d(da,f,u,appctx);CHKERRQ(ierr); // FIXME: Why does this give nonzero? 716*c4762a1bSJed Brown // } 717*c4762a1bSJed Brown 718*c4762a1bSJed Brown /* 719*c4762a1bSJed Brown Gather global vector, using the 2-step process 720*c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd(). 721*c4762a1bSJed Brown 722*c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or 723*c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above). 724*c4762a1bSJed Brown */ 725*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(da,localF,ADD_VALUES,F);CHKERRQ(ierr); 726*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(da,localF,ADD_VALUES,F);CHKERRQ(ierr); 727*c4762a1bSJed Brown 728*c4762a1bSJed Brown /* 729*c4762a1bSJed Brown Restore vectors 730*c4762a1bSJed Brown */ 731*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,localF,&f);CHKERRQ(ierr); 732*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 733*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr); 734*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 735*c4762a1bSJed Brown ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr); 736*c4762a1bSJed Brown 737*c4762a1bSJed Brown /* Destroy AFields appropriately */ 738*c4762a1bSJed Brown f_a += gys; 739*c4762a1bSJed Brown u_a += gys; 740*c4762a1bSJed Brown delete[] f_a; 741*c4762a1bSJed Brown delete[] u_a; 742*c4762a1bSJed Brown delete[] f_c; 743*c4762a1bSJed Brown delete[] u_c; 744*c4762a1bSJed Brown 745*c4762a1bSJed Brown PetscFunctionReturn(0); 746*c4762a1bSJed Brown } 747*c4762a1bSJed Brown 748*c4762a1bSJed Brown PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 749*c4762a1bSJed Brown { 750*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 751*c4762a1bSJed Brown DM da; 752*c4762a1bSJed Brown PetscErrorCode ierr; 753*c4762a1bSJed Brown PetscScalar *u_vec; 754*c4762a1bSJed Brown Vec localU; 755*c4762a1bSJed Brown 756*c4762a1bSJed Brown PetscFunctionBegin; 757*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 758*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 759*c4762a1bSJed Brown 760*c4762a1bSJed Brown /* 761*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 762*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 763*c4762a1bSJed Brown By placing code between these two statements, computations can be 764*c4762a1bSJed Brown done while messages are in transition. 765*c4762a1bSJed Brown */ 766*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 767*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 768*c4762a1bSJed Brown 769*c4762a1bSJed Brown /* Get pointers to vector data */ 770*c4762a1bSJed Brown ierr = VecGetArray(localU,&u_vec);CHKERRQ(ierr); 771*c4762a1bSJed Brown 772*c4762a1bSJed Brown /* 773*c4762a1bSJed Brown Compute Jacobian 774*c4762a1bSJed Brown */ 775*c4762a1bSJed Brown ierr = PetscAdolcComputeIJacobianLocalIDMass(1,A,u_vec,a,appctx->adctx);CHKERRQ(ierr); 776*c4762a1bSJed Brown 777*c4762a1bSJed Brown /* 778*c4762a1bSJed Brown Restore vectors 779*c4762a1bSJed Brown */ 780*c4762a1bSJed Brown ierr = VecRestoreArray(localU,&u_vec);CHKERRQ(ierr); 781*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 782*c4762a1bSJed Brown PetscFunctionReturn(0); 783*c4762a1bSJed Brown } 784*c4762a1bSJed Brown 785*c4762a1bSJed Brown PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 786*c4762a1bSJed Brown { 787*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 788*c4762a1bSJed Brown DM da; 789*c4762a1bSJed Brown PetscErrorCode ierr; 790*c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 791*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 792*c4762a1bSJed Brown PetscScalar uc,vc; 793*c4762a1bSJed Brown Field **u; 794*c4762a1bSJed Brown Vec localU; 795*c4762a1bSJed Brown MatStencil stencil[6],rowstencil; 796*c4762a1bSJed Brown PetscScalar entries[6]; 797*c4762a1bSJed Brown 798*c4762a1bSJed Brown PetscFunctionBegin; 799*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 800*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 801*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); 802*c4762a1bSJed Brown 803*c4762a1bSJed Brown hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 804*c4762a1bSJed Brown hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 805*c4762a1bSJed Brown 806*c4762a1bSJed Brown /* 807*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 808*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 809*c4762a1bSJed Brown By placing code between these two statements, computations can be 810*c4762a1bSJed Brown done while messages are in transition. 811*c4762a1bSJed Brown */ 812*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 813*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 814*c4762a1bSJed Brown 815*c4762a1bSJed Brown /* 816*c4762a1bSJed Brown Get pointers to vector data 817*c4762a1bSJed Brown */ 818*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 819*c4762a1bSJed Brown 820*c4762a1bSJed Brown /* 821*c4762a1bSJed Brown Get local grid boundaries 822*c4762a1bSJed Brown */ 823*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 824*c4762a1bSJed Brown 825*c4762a1bSJed Brown stencil[0].k = 0; 826*c4762a1bSJed Brown stencil[1].k = 0; 827*c4762a1bSJed Brown stencil[2].k = 0; 828*c4762a1bSJed Brown stencil[3].k = 0; 829*c4762a1bSJed Brown stencil[4].k = 0; 830*c4762a1bSJed Brown stencil[5].k = 0; 831*c4762a1bSJed Brown rowstencil.k = 0; 832*c4762a1bSJed Brown /* 833*c4762a1bSJed Brown Compute function over the locally owned part of the grid 834*c4762a1bSJed Brown */ 835*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 836*c4762a1bSJed Brown 837*c4762a1bSJed Brown stencil[0].j = j-1; 838*c4762a1bSJed Brown stencil[1].j = j+1; 839*c4762a1bSJed Brown stencil[2].j = j; 840*c4762a1bSJed Brown stencil[3].j = j; 841*c4762a1bSJed Brown stencil[4].j = j; 842*c4762a1bSJed Brown stencil[5].j = j; 843*c4762a1bSJed Brown rowstencil.k = 0; rowstencil.j = j; 844*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 845*c4762a1bSJed Brown uc = u[j][i].u; 846*c4762a1bSJed Brown vc = u[j][i].v; 847*c4762a1bSJed Brown 848*c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 849*c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 850*c4762a1bSJed Brown 851*c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 852*c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 853*c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 854*c4762a1bSJed Brown 855*c4762a1bSJed Brown stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy; 856*c4762a1bSJed Brown stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy; 857*c4762a1bSJed Brown stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx; 858*c4762a1bSJed Brown stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx; 859*c4762a1bSJed Brown stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a; 860*c4762a1bSJed Brown stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc; 861*c4762a1bSJed Brown rowstencil.i = i; rowstencil.c = 0; 862*c4762a1bSJed Brown 863*c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 864*c4762a1bSJed Brown if (appctx->aijpc) { 865*c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 866*c4762a1bSJed Brown } 867*c4762a1bSJed Brown stencil[0].c = 1; entries[0] = -appctx->D2*sy; 868*c4762a1bSJed Brown stencil[1].c = 1; entries[1] = -appctx->D2*sy; 869*c4762a1bSJed Brown stencil[2].c = 1; entries[2] = -appctx->D2*sx; 870*c4762a1bSJed Brown stencil[3].c = 1; entries[3] = -appctx->D2*sx; 871*c4762a1bSJed Brown stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a; 872*c4762a1bSJed Brown stencil[5].c = 0; entries[5] = -vc*vc; 873*c4762a1bSJed Brown 874*c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 875*c4762a1bSJed Brown if (appctx->aijpc) { 876*c4762a1bSJed Brown ierr = MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 877*c4762a1bSJed Brown } 878*c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 879*c4762a1bSJed Brown } 880*c4762a1bSJed Brown } 881*c4762a1bSJed Brown 882*c4762a1bSJed Brown /* 883*c4762a1bSJed Brown Restore vectors 884*c4762a1bSJed Brown */ 885*c4762a1bSJed Brown ierr = PetscLogFlops(19*xm*ym);CHKERRQ(ierr); 886*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 887*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 888*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 889*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 890*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 891*c4762a1bSJed Brown if (appctx->aijpc) { 892*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 893*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 894*c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 895*c4762a1bSJed Brown } 896*c4762a1bSJed Brown PetscFunctionReturn(0); 897*c4762a1bSJed Brown } 898*c4762a1bSJed Brown 899*c4762a1bSJed Brown PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 900*c4762a1bSJed Brown { 901*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 902*c4762a1bSJed Brown DM da; 903*c4762a1bSJed Brown PetscErrorCode ierr; 904*c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 905*c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 906*c4762a1bSJed Brown PetscScalar uc,vc; 907*c4762a1bSJed Brown Field **u; 908*c4762a1bSJed Brown Vec localU; 909*c4762a1bSJed Brown MatStencil stencil[6],rowstencil; 910*c4762a1bSJed Brown PetscScalar entries[6]; 911*c4762a1bSJed Brown 912*c4762a1bSJed Brown PetscFunctionBegin; 913*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 914*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 915*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); 916*c4762a1bSJed Brown 917*c4762a1bSJed Brown hx = 2.50/(PetscReal)(Mx); sx = 1.0/(hx*hx); 918*c4762a1bSJed Brown hy = 2.50/(PetscReal)(My); sy = 1.0/(hy*hy); 919*c4762a1bSJed Brown 920*c4762a1bSJed Brown /* 921*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 922*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 923*c4762a1bSJed Brown By placing code between these two statements, computations can be 924*c4762a1bSJed Brown done while messages are in transition. 925*c4762a1bSJed Brown */ 926*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 927*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 928*c4762a1bSJed Brown 929*c4762a1bSJed Brown /* 930*c4762a1bSJed Brown Get pointers to vector data 931*c4762a1bSJed Brown */ 932*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 933*c4762a1bSJed Brown 934*c4762a1bSJed Brown /* 935*c4762a1bSJed Brown Get local grid boundaries 936*c4762a1bSJed Brown */ 937*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 938*c4762a1bSJed Brown 939*c4762a1bSJed Brown stencil[0].k = 0; 940*c4762a1bSJed Brown stencil[1].k = 0; 941*c4762a1bSJed Brown stencil[2].k = 0; 942*c4762a1bSJed Brown stencil[3].k = 0; 943*c4762a1bSJed Brown stencil[4].k = 0; 944*c4762a1bSJed Brown stencil[5].k = 0; 945*c4762a1bSJed Brown rowstencil.k = 0; 946*c4762a1bSJed Brown 947*c4762a1bSJed Brown /* 948*c4762a1bSJed Brown Compute function over the locally owned part of the grid 949*c4762a1bSJed Brown */ 950*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 951*c4762a1bSJed Brown 952*c4762a1bSJed Brown stencil[0].j = j-1; 953*c4762a1bSJed Brown stencil[1].j = j+1; 954*c4762a1bSJed Brown stencil[2].j = j; 955*c4762a1bSJed Brown stencil[3].j = j; 956*c4762a1bSJed Brown stencil[4].j = j; 957*c4762a1bSJed Brown stencil[5].j = j; 958*c4762a1bSJed Brown rowstencil.k = 0; rowstencil.j = j; 959*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 960*c4762a1bSJed Brown uc = u[j][i].u; 961*c4762a1bSJed Brown vc = u[j][i].v; 962*c4762a1bSJed Brown 963*c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 964*c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 965*c4762a1bSJed Brown 966*c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 967*c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 968*c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 969*c4762a1bSJed Brown 970*c4762a1bSJed Brown stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy; 971*c4762a1bSJed Brown stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy; 972*c4762a1bSJed Brown stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx; 973*c4762a1bSJed Brown stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx; 974*c4762a1bSJed Brown stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma; 975*c4762a1bSJed Brown stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc; 976*c4762a1bSJed Brown rowstencil.i = i; rowstencil.c = 0; 977*c4762a1bSJed Brown 978*c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 979*c4762a1bSJed Brown 980*c4762a1bSJed Brown stencil[0].c = 1; entries[0] = appctx->D2*sy; 981*c4762a1bSJed Brown stencil[1].c = 1; entries[1] = appctx->D2*sy; 982*c4762a1bSJed Brown stencil[2].c = 1; entries[2] = appctx->D2*sx; 983*c4762a1bSJed Brown stencil[3].c = 1; entries[3] = appctx->D2*sx; 984*c4762a1bSJed Brown stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa; 985*c4762a1bSJed Brown stencil[5].c = 0; entries[5] = vc*vc; 986*c4762a1bSJed Brown rowstencil.c = 1; 987*c4762a1bSJed Brown 988*c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 989*c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 990*c4762a1bSJed Brown } 991*c4762a1bSJed Brown } 992*c4762a1bSJed Brown 993*c4762a1bSJed Brown /* 994*c4762a1bSJed Brown Restore vectors 995*c4762a1bSJed Brown */ 996*c4762a1bSJed Brown ierr = PetscLogFlops(19*xm*ym);CHKERRQ(ierr); 997*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 998*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 999*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1000*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1001*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1002*c4762a1bSJed Brown if (appctx->aijpc) { 1003*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1004*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1005*c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1006*c4762a1bSJed Brown } 1007*c4762a1bSJed Brown PetscFunctionReturn(0); 1008*c4762a1bSJed Brown } 1009*c4762a1bSJed Brown 1010*c4762a1bSJed Brown PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 1011*c4762a1bSJed Brown { 1012*c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 1013*c4762a1bSJed Brown DM da; 1014*c4762a1bSJed Brown PetscErrorCode ierr; 1015*c4762a1bSJed Brown PetscScalar *u_vec; 1016*c4762a1bSJed Brown Vec localU; 1017*c4762a1bSJed Brown 1018*c4762a1bSJed Brown PetscFunctionBegin; 1019*c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 1020*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 1021*c4762a1bSJed Brown 1022*c4762a1bSJed Brown /* 1023*c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 1024*c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 1025*c4762a1bSJed Brown By placing code between these two statements, computations can be 1026*c4762a1bSJed Brown done while messages are in transition. 1027*c4762a1bSJed Brown */ 1028*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 1029*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 1030*c4762a1bSJed Brown 1031*c4762a1bSJed Brown /* Get pointers to vector data */ 1032*c4762a1bSJed Brown ierr = VecGetArray(localU,&u_vec);CHKERRQ(ierr); 1033*c4762a1bSJed Brown 1034*c4762a1bSJed Brown /* 1035*c4762a1bSJed Brown Compute Jacobian 1036*c4762a1bSJed Brown */ 1037*c4762a1bSJed Brown ierr = PetscAdolcComputeRHSJacobianLocal(1,A,u_vec,appctx->adctx);CHKERRQ(ierr); 1038*c4762a1bSJed Brown 1039*c4762a1bSJed Brown /* 1040*c4762a1bSJed Brown Restore vectors 1041*c4762a1bSJed Brown */ 1042*c4762a1bSJed Brown ierr = VecRestoreArray(localU,&u_vec);CHKERRQ(ierr); 1043*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 1044*c4762a1bSJed Brown PetscFunctionReturn(0); 1045*c4762a1bSJed Brown } 1046*c4762a1bSJed Brown 1047*c4762a1bSJed Brown /*TEST 1048*c4762a1bSJed Brown 1049*c4762a1bSJed Brown build: 1050*c4762a1bSJed Brown requires: double !complex adolc colpack 1051*c4762a1bSJed Brown 1052*c4762a1bSJed Brown test: 1053*c4762a1bSJed Brown suffix: 1 1054*c4762a1bSJed Brown nsize: 1 1055*c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian 1056*c4762a1bSJed Brown output_file: output/adr_ex5adj_1.out 1057*c4762a1bSJed Brown 1058*c4762a1bSJed Brown test: 1059*c4762a1bSJed Brown suffix: 2 1060*c4762a1bSJed Brown nsize: 1 1061*c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian -implicitform 1062*c4762a1bSJed Brown output_file: output/adr_ex5adj_2.out 1063*c4762a1bSJed Brown 1064*c4762a1bSJed Brown test: 1065*c4762a1bSJed Brown suffix: 3 1066*c4762a1bSJed Brown nsize: 4 1067*c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor 1068*c4762a1bSJed Brown output_file: output/adr_ex5adj_3.out 1069*c4762a1bSJed Brown 1070*c4762a1bSJed Brown test: 1071*c4762a1bSJed Brown suffix: 4 1072*c4762a1bSJed Brown nsize: 4 1073*c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor -implicitform 1074*c4762a1bSJed Brown output_file: output/adr_ex5adj_4.out 1075*c4762a1bSJed Brown 1076*c4762a1bSJed Brown testset: 1077*c4762a1bSJed Brown suffix: 5 1078*c4762a1bSJed Brown nsize: 4 1079*c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse 1080*c4762a1bSJed Brown output_file: output/adr_ex5adj_5.out 1081*c4762a1bSJed Brown 1082*c4762a1bSJed Brown testset: 1083*c4762a1bSJed Brown suffix: 6 1084*c4762a1bSJed Brown nsize: 4 1085*c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse -implicitform 1086*c4762a1bSJed Brown output_file: output/adr_ex5adj_6.out 1087*c4762a1bSJed Brown 1088*c4762a1bSJed Brown TEST*/ 1089