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