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