1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown See ex5.c for details on the equation. 5c4762a1bSJed Brown This code demonestrates the TSAdjoint interface to a system of time-dependent partial differential equations. 6c4762a1bSJed Brown It computes the sensitivity of a component in the final solution, which locates in the center of the 2D domain, w.r.t. the initial conditions. 7c4762a1bSJed Brown The user does not need to provide any additional functions. The required functions in the original simulation are reused in the adjoint run. 8c4762a1bSJed Brown 9c4762a1bSJed Brown Runtime options: 10c4762a1bSJed Brown -forwardonly - run the forward simulation without adjoint 11c4762a1bSJed Brown -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used 12c4762a1bSJed Brown -aijpc - set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL) 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscsys.h> 16c4762a1bSJed Brown #include <petscdm.h> 17c4762a1bSJed Brown #include <petscdmda.h> 18c4762a1bSJed Brown #include <petscts.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown typedef struct { 21c4762a1bSJed Brown PetscScalar u,v; 22c4762a1bSJed Brown } Field; 23c4762a1bSJed Brown 24c4762a1bSJed Brown typedef struct { 25c4762a1bSJed Brown PetscReal D1,D2,gamma,kappa; 26c4762a1bSJed Brown PetscBool aijpc; 27c4762a1bSJed Brown } AppCtx; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* 30c4762a1bSJed Brown User-defined routines 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 33c4762a1bSJed Brown PetscErrorCode InitialConditions(DM,Vec); 34c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 35c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 36c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 41c4762a1bSJed Brown PetscErrorCode ierr; 42c4762a1bSJed Brown Field **l; 43c4762a1bSJed Brown PetscFunctionBegin; 44c4762a1bSJed Brown 45c4762a1bSJed 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); 46c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 47c4762a1bSJed Brown i = (PetscInt)(x*(Mx-1)); 48c4762a1bSJed Brown j = (PetscInt)(y*(My-1)); 49c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) { 52c4762a1bSJed Brown /* the i,j vertex is on this process */ 53c4762a1bSJed Brown ierr = DMDAVecGetArray(da,lambda,&l);CHKERRQ(ierr); 54c4762a1bSJed Brown l[j][i].u = 1.0; 55c4762a1bSJed Brown l[j][i].v = 1.0; 56c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,lambda,&l);CHKERRQ(ierr); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown int main(int argc,char **argv) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown TS ts; /* ODE integrator */ 64c4762a1bSJed Brown Vec x; /* solution */ 65c4762a1bSJed Brown PetscErrorCode ierr; 66c4762a1bSJed Brown DM da; 67c4762a1bSJed Brown AppCtx appctx; 68c4762a1bSJed Brown Vec lambda[1]; 69c4762a1bSJed Brown PetscBool forwardonly=PETSC_FALSE,implicitform=PETSC_TRUE; 70c4762a1bSJed Brown 71c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 72c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr); 74c4762a1bSJed Brown appctx.aijpc = PETSC_FALSE; 75c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown appctx.D1 = 8.0e-5; 78c4762a1bSJed Brown appctx.D2 = 4.0e-5; 79c4762a1bSJed Brown appctx.gamma = .024; 80c4762a1bSJed Brown appctx.kappa = .06; 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 85c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 93c4762a1bSJed Brown vectors that are the same types 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 95c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98c4762a1bSJed Brown Create timestepping solver context 99c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 100c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 104c4762a1bSJed Brown if (!implicitform) { 105c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr); 108c4762a1bSJed Brown } else { 109c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&appctx);CHKERRQ(ierr); 111c4762a1bSJed Brown if (appctx.aijpc) { 112c4762a1bSJed Brown Mat A,B; 113c4762a1bSJed Brown 114c4762a1bSJed Brown ierr = DMSetMatType(da,MATSELL);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 117c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */ 118c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,B,IJacobian,&appctx);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 121c4762a1bSJed Brown } else { 122c4762a1bSJed Brown ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,&appctx);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Set initial conditions 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129c4762a1bSJed Brown ierr = InitialConditions(da,x);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* 133c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 134c4762a1bSJed Brown */ 135c4762a1bSJed Brown if (!forwardonly) { ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Set solver options 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 140c4762a1bSJed Brown ierr = TSSetMaxTime(ts,200.0);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.5);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Solve ODE system 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 149c4762a1bSJed Brown if (!forwardonly) { 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Start the Adjoint model 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153c4762a1bSJed Brown ierr = VecDuplicate(x,&lambda[0]);CHKERRQ(ierr); 154c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 155c4762a1bSJed Brown ierr = InitializeLambda(da,lambda[0],0.5,0.5);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,lambda,NULL);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 162c4762a1bSJed Brown are no longer needed. 163c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 164c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = PetscFinalize(); 168c4762a1bSJed Brown return ierr; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown RHSFunction - Evaluates nonlinear function, F(x). 174c4762a1bSJed Brown 175c4762a1bSJed Brown Input Parameters: 176c4762a1bSJed Brown . ts - the TS context 177c4762a1bSJed Brown . X - input vector 178c4762a1bSJed Brown . ptr - optional user-defined context, as set by TSSetRHSFunction() 179c4762a1bSJed Brown 180c4762a1bSJed Brown Output Parameter: 181c4762a1bSJed Brown . F - function vector 182c4762a1bSJed Brown */ 183c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 186c4762a1bSJed Brown DM da; 187c4762a1bSJed Brown PetscErrorCode ierr; 188c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 189c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 190c4762a1bSJed Brown PetscScalar uc,uxx,uyy,vc,vxx,vyy; 191c4762a1bSJed Brown Field **u,**f; 192c4762a1bSJed Brown Vec localU; 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscFunctionBegin; 195c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 197c4762a1bSJed 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); 198c4762a1bSJed Brown hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 199c4762a1bSJed Brown hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 202c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 203c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 204c4762a1bSJed Brown By placing code between these two statements, computations can be 205c4762a1bSJed Brown done while messages are in transition. 206c4762a1bSJed Brown */ 207c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Get pointers to vector data 212c4762a1bSJed Brown */ 213c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* 217c4762a1bSJed Brown Get local grid boundaries 218c4762a1bSJed Brown */ 219c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* 222c4762a1bSJed Brown Compute function over the locally owned part of the grid 223c4762a1bSJed Brown */ 224c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 225c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 226c4762a1bSJed Brown uc = u[j][i].u; 227c4762a1bSJed Brown uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 228c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 229c4762a1bSJed Brown vc = u[j][i].v; 230c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 231c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 232c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc); 233c4762a1bSJed Brown f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 236ca0c957dSBarry Smith ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Restore vectors 240c4762a1bSJed Brown */ 241c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 244c4762a1bSJed Brown PetscFunctionReturn(0); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 248c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown PetscErrorCode ierr; 251c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym,Mx,My; 252c4762a1bSJed Brown Field **u; 253c4762a1bSJed Brown PetscReal hx,hy,x,y; 254c4762a1bSJed Brown 255c4762a1bSJed Brown PetscFunctionBegin; 256c4762a1bSJed 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); 257c4762a1bSJed Brown 258c4762a1bSJed Brown hx = 2.5/(PetscReal)Mx; 259c4762a1bSJed Brown hy = 2.5/(PetscReal)My; 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* 262c4762a1bSJed Brown Get pointers to vector data 263c4762a1bSJed Brown */ 264c4762a1bSJed Brown ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* 267c4762a1bSJed Brown Get local grid boundaries 268c4762a1bSJed Brown */ 269c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* 272c4762a1bSJed Brown Compute function over the locally owned part of the grid 273c4762a1bSJed Brown */ 274c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 275c4762a1bSJed Brown y = j*hy; 276c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 277c4762a1bSJed Brown x = i*hx; 278c4762a1bSJed 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); 279c4762a1bSJed Brown else u[j][i].v = 0.0; 280c4762a1bSJed Brown 281c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0*u[j][i].v; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* 286c4762a1bSJed Brown Restore vectors 287c4762a1bSJed Brown */ 288c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 289c4762a1bSJed Brown PetscFunctionReturn(0); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx) 293c4762a1bSJed Brown { 294c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 295c4762a1bSJed Brown DM da; 296c4762a1bSJed Brown PetscErrorCode ierr; 297c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 298c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 299c4762a1bSJed Brown PetscScalar uc,vc; 300c4762a1bSJed Brown Field **u; 301c4762a1bSJed Brown Vec localU; 302c4762a1bSJed Brown MatStencil stencil[6],rowstencil; 303c4762a1bSJed Brown PetscScalar entries[6]; 304c4762a1bSJed Brown 305c4762a1bSJed Brown PetscFunctionBegin; 306c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 308c4762a1bSJed 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); 309c4762a1bSJed Brown 310c4762a1bSJed Brown hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 311c4762a1bSJed Brown hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 315c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 316c4762a1bSJed Brown By placing code between these two statements, computations can be 317c4762a1bSJed Brown done while messages are in transition. 318c4762a1bSJed Brown */ 319c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 321c4762a1bSJed Brown 322c4762a1bSJed Brown /* 323c4762a1bSJed Brown Get pointers to vector data 324c4762a1bSJed Brown */ 325c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* 328c4762a1bSJed Brown Get local grid boundaries 329c4762a1bSJed Brown */ 330c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 331c4762a1bSJed Brown 332c4762a1bSJed Brown stencil[0].k = 0; 333c4762a1bSJed Brown stencil[1].k = 0; 334c4762a1bSJed Brown stencil[2].k = 0; 335c4762a1bSJed Brown stencil[3].k = 0; 336c4762a1bSJed Brown stencil[4].k = 0; 337c4762a1bSJed Brown stencil[5].k = 0; 338c4762a1bSJed Brown rowstencil.k = 0; 339c4762a1bSJed Brown /* 340c4762a1bSJed Brown Compute function over the locally owned part of the grid 341c4762a1bSJed Brown */ 342c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 343c4762a1bSJed Brown 344c4762a1bSJed Brown stencil[0].j = j-1; 345c4762a1bSJed Brown stencil[1].j = j+1; 346c4762a1bSJed Brown stencil[2].j = j; 347c4762a1bSJed Brown stencil[3].j = j; 348c4762a1bSJed Brown stencil[4].j = j; 349c4762a1bSJed Brown stencil[5].j = j; 350c4762a1bSJed Brown rowstencil.k = 0; rowstencil.j = j; 351c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 352c4762a1bSJed Brown uc = u[j][i].u; 353c4762a1bSJed Brown vc = u[j][i].v; 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 356c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 357c4762a1bSJed Brown 358c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 359c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 360c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 361c4762a1bSJed Brown 362c4762a1bSJed Brown stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy; 363c4762a1bSJed Brown stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy; 364c4762a1bSJed Brown stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx; 365c4762a1bSJed Brown stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx; 366c4762a1bSJed Brown stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma; 367c4762a1bSJed Brown stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc; 368c4762a1bSJed Brown rowstencil.i = i; rowstencil.c = 0; 369c4762a1bSJed Brown 370c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 371c4762a1bSJed Brown if (appctx->aijpc) { 372c4762a1bSJed Brown ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown stencil[0].c = 1; entries[0] = appctx->D2*sy; 375c4762a1bSJed Brown stencil[1].c = 1; entries[1] = appctx->D2*sy; 376c4762a1bSJed Brown stencil[2].c = 1; entries[2] = appctx->D2*sx; 377c4762a1bSJed Brown stencil[3].c = 1; entries[3] = appctx->D2*sx; 378c4762a1bSJed Brown stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa; 379c4762a1bSJed Brown stencil[5].c = 0; entries[5] = vc*vc; 380c4762a1bSJed Brown rowstencil.c = 1; 381c4762a1bSJed Brown 382c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 383c4762a1bSJed Brown if (appctx->aijpc) { 384c4762a1bSJed Brown ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 387c4762a1bSJed Brown } 388c4762a1bSJed Brown } 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* 391c4762a1bSJed Brown Restore vectors 392c4762a1bSJed Brown */ 393ca0c957dSBarry Smith ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 399c4762a1bSJed Brown if (appctx->aijpc) { 400c4762a1bSJed Brown ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 401c4762a1bSJed Brown ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 402c4762a1bSJed Brown ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown PetscFunctionReturn(0); 406c4762a1bSJed Brown } 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* 409c4762a1bSJed Brown IFunction - Evaluates implicit nonlinear function, xdot - F(x). 410c4762a1bSJed Brown 411c4762a1bSJed Brown Input Parameters: 412c4762a1bSJed Brown . ts - the TS context 413c4762a1bSJed Brown . U - input vector 414c4762a1bSJed Brown . Udot - input vector 415c4762a1bSJed Brown . ptr - optional user-defined context, as set by TSSetRHSFunction() 416c4762a1bSJed Brown 417c4762a1bSJed Brown Output Parameter: 418c4762a1bSJed Brown . F - function vector 419c4762a1bSJed Brown */ 420c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr) 421c4762a1bSJed Brown { 422c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ptr; 423c4762a1bSJed Brown DM da; 424c4762a1bSJed Brown PetscErrorCode ierr; 425c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 426c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 427c4762a1bSJed Brown PetscScalar uc,uxx,uyy,vc,vxx,vyy; 428c4762a1bSJed Brown Field **u,**f,**udot; 429c4762a1bSJed Brown Vec localU; 430c4762a1bSJed Brown 431c4762a1bSJed Brown PetscFunctionBegin; 432c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 433c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 434c4762a1bSJed 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); 435c4762a1bSJed Brown hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 436c4762a1bSJed Brown hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 437c4762a1bSJed Brown 438c4762a1bSJed Brown /* 439c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 440c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 441c4762a1bSJed Brown By placing code between these two statements, computations can be 442c4762a1bSJed Brown done while messages are in transition. 443c4762a1bSJed Brown */ 444c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 446c4762a1bSJed Brown 447c4762a1bSJed Brown /* 448c4762a1bSJed Brown Get pointers to vector data 449c4762a1bSJed Brown */ 450c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 451c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr); 452c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr); 453c4762a1bSJed Brown 454c4762a1bSJed Brown /* 455c4762a1bSJed Brown Get local grid boundaries 456c4762a1bSJed Brown */ 457c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 458c4762a1bSJed Brown 459c4762a1bSJed Brown /* 460c4762a1bSJed Brown Compute function over the locally owned part of the grid 461c4762a1bSJed Brown */ 462c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 463c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 464c4762a1bSJed Brown uc = u[j][i].u; 465c4762a1bSJed Brown uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 466c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 467c4762a1bSJed Brown vc = u[j][i].v; 468c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 469c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 470c4762a1bSJed Brown f[j][i].u = udot[j][i].u - ( appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc)); 471c4762a1bSJed Brown f[j][i].v = udot[j][i].v - ( appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc); 472c4762a1bSJed Brown } 473c4762a1bSJed Brown } 474ca0c957dSBarry Smith ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr); 475c4762a1bSJed Brown 476c4762a1bSJed Brown /* 477c4762a1bSJed Brown Restore vectors 478c4762a1bSJed Brown */ 479c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 480c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr); 481c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 483c4762a1bSJed Brown PetscFunctionReturn(0); 484c4762a1bSJed Brown } 485c4762a1bSJed Brown 486c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx) 487c4762a1bSJed Brown { 488c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 489c4762a1bSJed Brown DM da; 490c4762a1bSJed Brown PetscErrorCode ierr; 491c4762a1bSJed Brown PetscInt i,j,Mx,My,xs,ys,xm,ym; 492c4762a1bSJed Brown PetscReal hx,hy,sx,sy; 493c4762a1bSJed Brown PetscScalar uc,vc; 494c4762a1bSJed Brown Field **u; 495c4762a1bSJed Brown Vec localU; 496c4762a1bSJed Brown MatStencil stencil[6],rowstencil; 497c4762a1bSJed Brown PetscScalar entries[6]; 498c4762a1bSJed Brown 499c4762a1bSJed Brown PetscFunctionBegin; 500c4762a1bSJed Brown ierr = TSGetDM(ts,&da);CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr); 502c4762a1bSJed 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); 503c4762a1bSJed Brown 504c4762a1bSJed Brown hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 505c4762a1bSJed Brown hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 506c4762a1bSJed Brown 507c4762a1bSJed Brown /* 508c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 509c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 510c4762a1bSJed Brown By placing code between these two statements, computations can be 511c4762a1bSJed Brown done while messages are in transition. 512c4762a1bSJed Brown */ 513c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 514c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr); 515c4762a1bSJed Brown 516c4762a1bSJed Brown /* 517c4762a1bSJed Brown Get pointers to vector data 518c4762a1bSJed Brown */ 519c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr); 520c4762a1bSJed Brown 521c4762a1bSJed Brown /* 522c4762a1bSJed Brown Get local grid boundaries 523c4762a1bSJed Brown */ 524c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 525c4762a1bSJed Brown 526c4762a1bSJed Brown stencil[0].k = 0; 527c4762a1bSJed Brown stencil[1].k = 0; 528c4762a1bSJed Brown stencil[2].k = 0; 529c4762a1bSJed Brown stencil[3].k = 0; 530c4762a1bSJed Brown stencil[4].k = 0; 531c4762a1bSJed Brown stencil[5].k = 0; 532c4762a1bSJed Brown rowstencil.k = 0; 533c4762a1bSJed Brown /* 534c4762a1bSJed Brown Compute function over the locally owned part of the grid 535c4762a1bSJed Brown */ 536c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 537c4762a1bSJed Brown 538c4762a1bSJed Brown stencil[0].j = j-1; 539c4762a1bSJed Brown stencil[1].j = j+1; 540c4762a1bSJed Brown stencil[2].j = j; 541c4762a1bSJed Brown stencil[3].j = j; 542c4762a1bSJed Brown stencil[4].j = j; 543c4762a1bSJed Brown stencil[5].j = j; 544c4762a1bSJed Brown rowstencil.k = 0; rowstencil.j = j; 545c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 546c4762a1bSJed Brown uc = u[j][i].u; 547c4762a1bSJed Brown vc = u[j][i].v; 548c4762a1bSJed Brown 549c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 550c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 551c4762a1bSJed Brown 552c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 553c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 554c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 555c4762a1bSJed Brown 556c4762a1bSJed Brown stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy; 557c4762a1bSJed Brown stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy; 558c4762a1bSJed Brown stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx; 559c4762a1bSJed Brown stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx; 560c4762a1bSJed Brown stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a; 561c4762a1bSJed Brown stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc; 562c4762a1bSJed Brown rowstencil.i = i; rowstencil.c = 0; 563c4762a1bSJed Brown 564c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 565c4762a1bSJed Brown if (appctx->aijpc) { 566c4762a1bSJed Brown ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 567c4762a1bSJed Brown } 568c4762a1bSJed Brown stencil[0].c = 1; entries[0] = -appctx->D2*sy; 569c4762a1bSJed Brown stencil[1].c = 1; entries[1] = -appctx->D2*sy; 570c4762a1bSJed Brown stencil[2].c = 1; entries[2] = -appctx->D2*sx; 571c4762a1bSJed Brown stencil[3].c = 1; entries[3] = -appctx->D2*sx; 572c4762a1bSJed Brown stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a; 573c4762a1bSJed Brown stencil[5].c = 0; entries[5] = -vc*vc; 574c4762a1bSJed Brown rowstencil.c = 1; 575c4762a1bSJed Brown 576c4762a1bSJed Brown ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 577c4762a1bSJed Brown if (appctx->aijpc) { 578c4762a1bSJed Brown ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } 583c4762a1bSJed Brown 584c4762a1bSJed Brown /* 585c4762a1bSJed Brown Restore vectors 586c4762a1bSJed Brown */ 587ca0c957dSBarry Smith ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr); 588c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr); 589c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr); 590c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 591c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 593c4762a1bSJed Brown if (appctx->aijpc) { 594c4762a1bSJed Brown ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 595c4762a1bSJed Brown ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596c4762a1bSJed Brown ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 597c4762a1bSJed Brown } 598c4762a1bSJed Brown PetscFunctionReturn(0); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown 601c4762a1bSJed Brown 602c4762a1bSJed Brown /*TEST 603c4762a1bSJed Brown 604c4762a1bSJed Brown build: 605c4762a1bSJed Brown requires: !complex !single 606c4762a1bSJed Brown 607c4762a1bSJed Brown test: 608c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20 609c4762a1bSJed Brown output_file: output/ex5adj_1.out 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 2 613c4762a1bSJed Brown nsize: 2 614c4762a1bSJed Brown args: -ts_max_steps 10 -ts_dt 10 -ts_monitor -ts_adjoint_monitor -ksp_monitor_short -da_grid_x 20 -da_grid_y 20 -ts_trajectory_dirname Test-dir -ts_trajectory_file_template test-%06D.cp 615c4762a1bSJed Brown 616c4762a1bSJed Brown test: 617c4762a1bSJed Brown suffix: 3 618c4762a1bSJed Brown nsize: 2 619c4762a1bSJed Brown args: -ts_max_steps 10 -ts_dt 10 -ts_adjoint_monitor_draw_sensi 620c4762a1bSJed Brown 621c4762a1bSJed Brown test: 622c4762a1bSJed Brown suffix: 4 623c4762a1bSJed Brown nsize: 2 624c4762a1bSJed Brown args: -ts_max_steps 10 -ts_dt 10 -ts_monitor -ts_adjoint_monitor -ksp_monitor_short -da_grid_x 20 -da_grid_y 20 -snes_fd_color 625c4762a1bSJed Brown output_file: output/ex5adj_2.out 626c4762a1bSJed Brown 627c4762a1bSJed Brown test: 628c4762a1bSJed Brown suffix: 5 629c4762a1bSJed Brown nsize: 2 630c4762a1bSJed Brown args: -ts_max_steps 10 -implicitform 0 -ts_type rk -ts_rk_type 4 -ts_monitor -ts_adjoint_monitor -da_grid_x 20 -da_grid_y 20 -snes_fd_color 631c4762a1bSJed Brown output_file: output/ex5adj_1.out 632c4762a1bSJed Brown 633c4762a1bSJed Brown test: 634c4762a1bSJed Brown suffix: knl 635c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -malloc_hbw -ts_trajectory_use_dram 1 636c4762a1bSJed Brown output_file: output/ex5adj_3.out 637c4762a1bSJed Brown requires: knl 638c4762a1bSJed Brown 639c4762a1bSJed Brown test: 640c4762a1bSJed Brown suffix: sell 641c4762a1bSJed Brown nsize: 4 642c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type none 643c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out 644c4762a1bSJed Brown 645c4762a1bSJed Brown test: 646c4762a1bSJed Brown suffix: aijsell 647c4762a1bSJed Brown nsize: 4 648c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type none 649c4762a1bSJed Brown output_file: output/ex5adj_sell_1.out 650c4762a1bSJed Brown 651c4762a1bSJed Brown test: 652c4762a1bSJed Brown suffix: sell2 653c4762a1bSJed Brown nsize: 4 654c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type sor 655c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out 656c4762a1bSJed Brown 657c4762a1bSJed Brown test: 658c4762a1bSJed Brown suffix: aijsell2 659c4762a1bSJed Brown nsize: 4 660c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type sor 661c4762a1bSJed Brown output_file: output/ex5adj_sell_2.out 662c4762a1bSJed Brown 663c4762a1bSJed Brown test: 664c4762a1bSJed Brown suffix: sell3 665c4762a1bSJed Brown nsize: 4 666c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bjacobi -mg_levels_pc_type bjacobi 667c4762a1bSJed Brown output_file: output/ex5adj_sell_3.out 668c4762a1bSJed Brown 669c4762a1bSJed Brown test: 670c4762a1bSJed Brown suffix: sell4 671c4762a1bSJed Brown nsize: 4 672c4762a1bSJed Brown args: -forwardonly -implicitform -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bjacobi -mg_levels_pc_type bjacobi 673c4762a1bSJed Brown output_file: output/ex5adj_sell_4.out 674c4762a1bSJed Brown 675c4762a1bSJed Brown test: 676c4762a1bSJed Brown suffix: sell5 677c4762a1bSJed Brown nsize: 4 678c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type sell -aijpc 679c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out 680c4762a1bSJed Brown 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: aijsell5 683c4762a1bSJed Brown nsize: 4 684c4762a1bSJed Brown args: -forwardonly -ts_max_steps 10 -ts_monitor -snes_monitor_short -dm_mat_type aijsell 685c4762a1bSJed Brown output_file: output/ex5adj_sell_5.out 686c4762a1bSJed Brown 687c4762a1bSJed Brown test: 688c4762a1bSJed Brown suffix: sell6 689c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -dm_mat_type sell -pc_type jacobi 690c4762a1bSJed Brown output_file: output/ex5adj_sell_6.out 691c4762a1bSJed Brown 692*fcb023d4SJed Brown test: 693*fcb023d4SJed Brown suffix: gamg1 694*fcb023d4SJed Brown args: -pc_type gamg -pc_mg_levels 2 -mg_levels_pc_type jacobi -ksp_monitor_short -ts_max_steps 2 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory 695*fcb023d4SJed Brown output_file: output/ex5adj_gamg_1.out 696*fcb023d4SJed Brown 697c4762a1bSJed Brown TEST*/ 698