160f0b76eSHong Zhang static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n"; 260f0b76eSHong Zhang 360f0b76eSHong Zhang /* 460f0b76eSHong Zhang See ex5.c for details on the equation. 560f0b76eSHong Zhang This code demonestrates the TSAdjoint interface to a system of time-dependent partial differential equations. 660f0b76eSHong Zhang 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. 760f0b76eSHong Zhang The user does not need to provide any additional functions. The required functions in the original simulation are reused in the adjoint run. 860f0b76eSHong Zhang 960f0b76eSHong Zhang Runtime options: 1060f0b76eSHong Zhang -forwardonly - run the forward simulation without adjoint 1160f0b76eSHong Zhang -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used 1260f0b76eSHong Zhang -aijpc - set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL) 1360f0b76eSHong Zhang */ 1460f0b76eSHong Zhang 1560f0b76eSHong Zhang #include "reaction_diffusion.h" 1660f0b76eSHong Zhang #include <petscdm.h> 1760f0b76eSHong Zhang #include <petscdmda.h> 1860f0b76eSHong Zhang 1960f0b76eSHong Zhang /* Matrix free support */ 2060f0b76eSHong Zhang typedef struct { 2160f0b76eSHong Zhang PetscReal time; 2260f0b76eSHong Zhang Vec U; 2360f0b76eSHong Zhang Vec Udot; 2460f0b76eSHong Zhang PetscReal shift; 2560f0b76eSHong Zhang AppCtx* appctx; 2660f0b76eSHong Zhang TS ts; 2760f0b76eSHong Zhang } MCtx; 2860f0b76eSHong Zhang 2960f0b76eSHong Zhang /* 3060f0b76eSHong Zhang User-defined routines 3160f0b76eSHong Zhang */ 3260f0b76eSHong Zhang PetscErrorCode InitialConditions(DM,Vec); 3360f0b76eSHong Zhang PetscErrorCode RHSJacobianShell(TS,PetscReal,Vec,Mat,Mat,void*); 3460f0b76eSHong Zhang PetscErrorCode IJacobianShell(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 3560f0b76eSHong Zhang 3660f0b76eSHong Zhang PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y) 3760f0b76eSHong Zhang { 3860f0b76eSHong Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 3960f0b76eSHong Zhang Field **l; 4060f0b76eSHong Zhang PetscFunctionBegin; 4160f0b76eSHong Zhang 42*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)); 4360f0b76eSHong Zhang /* locate the global i index for x and j index for y */ 4460f0b76eSHong Zhang i = (PetscInt)(x*(Mx-1)); 4560f0b76eSHong Zhang j = (PetscInt)(y*(My-1)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 4760f0b76eSHong Zhang 4860f0b76eSHong Zhang if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) { 4960f0b76eSHong Zhang /* the i,j vertex is on this process */ 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,lambda,&l)); 5160f0b76eSHong Zhang l[j][i].u = 1.0; 5260f0b76eSHong Zhang l[j][i].v = 1.0; 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,lambda,&l)); 5460f0b76eSHong Zhang } 5560f0b76eSHong Zhang PetscFunctionReturn(0); 5660f0b76eSHong Zhang } 5760f0b76eSHong Zhang 5860f0b76eSHong Zhang static PetscErrorCode MyRHSMatMultTranspose(Mat A_shell,Vec X,Vec Y) 5960f0b76eSHong Zhang { 6060f0b76eSHong Zhang MCtx *mctx; 6160f0b76eSHong Zhang AppCtx *appctx; 6260f0b76eSHong Zhang DM da; 6360f0b76eSHong Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 6460f0b76eSHong Zhang PetscReal hx,hy,sx,sy; 6560f0b76eSHong Zhang PetscScalar uc,uxx,uyy,vc,vxx,vyy,ucb,vcb; 6660f0b76eSHong Zhang Field **u,**x,**y; 6760f0b76eSHong Zhang Vec localX; 6860f0b76eSHong Zhang 6960f0b76eSHong Zhang PetscFunctionBeginUser; 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A_shell,&mctx)); 7160f0b76eSHong Zhang appctx = mctx->appctx; 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(mctx->ts,&da)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 74*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)); 7560f0b76eSHong Zhang hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 7660f0b76eSHong Zhang hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 7760f0b76eSHong Zhang 7860f0b76eSHong Zhang /* Scatter ghost points to local vector,using the 2-step process 7960f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 8060f0b76eSHong Zhang By placing code between these two statements, computations can be 8160f0b76eSHong Zhang done while messages are in transition. */ 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 8460f0b76eSHong Zhang 8560f0b76eSHong Zhang /* Get pointers to vector data */ 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localX,&x)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,mctx->U,&u)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Y,&y)); 8960f0b76eSHong Zhang 9060f0b76eSHong Zhang /* Get local grid boundaries */ 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 9260f0b76eSHong Zhang 9360f0b76eSHong Zhang /* Compute function over the locally owned part of the grid */ 9460f0b76eSHong Zhang for (j=ys; j<ys+ym; j++) { 9560f0b76eSHong Zhang for (i=xs; i<xs+xm; i++) { 9660f0b76eSHong Zhang uc = u[j][i].u; 9760f0b76eSHong Zhang ucb = x[j][i].u; 9860f0b76eSHong Zhang uxx = (-2.0*ucb + x[j][i-1].u + x[j][i+1].u)*sx; 9960f0b76eSHong Zhang uyy = (-2.0*ucb + x[j-1][i].u + x[j+1][i].u)*sy; 10060f0b76eSHong Zhang vc = u[j][i].v; 10160f0b76eSHong Zhang vcb = x[j][i].v; 10260f0b76eSHong Zhang vxx = (-2.0*vcb + x[j][i-1].v + x[j][i+1].v)*sx; 10360f0b76eSHong Zhang vyy = (-2.0*vcb + x[j-1][i].v + x[j+1][i].v)*sy; 10460f0b76eSHong Zhang y[j][i].u = appctx->D1*(uxx + uyy) - ucb*(vc*vc+appctx->gamma) + vc*vc*vcb; 10560f0b76eSHong Zhang y[j][i].v = appctx->D2*(vxx + vyy) - 2.0*uc*vc*ucb + (2.0*uc*vc-appctx->gamma - appctx->kappa)*vcb; 10660f0b76eSHong Zhang } 10760f0b76eSHong Zhang } 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,mctx->U,&u)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Y,&y)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 11260f0b76eSHong Zhang PetscFunctionReturn(0); 11360f0b76eSHong Zhang } 11460f0b76eSHong Zhang 11560f0b76eSHong Zhang static PetscErrorCode MyIMatMultTranspose(Mat A_shell,Vec X,Vec Y) 11660f0b76eSHong Zhang { 11760f0b76eSHong Zhang MCtx *mctx; 11860f0b76eSHong Zhang AppCtx *appctx; 11960f0b76eSHong Zhang DM da; 12060f0b76eSHong Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 12160f0b76eSHong Zhang PetscReal hx,hy,sx,sy; 12260f0b76eSHong Zhang PetscScalar uc,uxx,uyy,vc,vxx,vyy,ucb,vcb; 12360f0b76eSHong Zhang Field **u,**x,**y; 12460f0b76eSHong Zhang Vec localX; 12560f0b76eSHong Zhang 12660f0b76eSHong Zhang PetscFunctionBeginUser; 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A_shell,&mctx)); 12860f0b76eSHong Zhang appctx = mctx->appctx; 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(mctx->ts,&da)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 131*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)); 13260f0b76eSHong Zhang hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 13360f0b76eSHong Zhang hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 13460f0b76eSHong Zhang 13560f0b76eSHong Zhang /* Scatter ghost points to local vector,using the 2-step process 13660f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 13760f0b76eSHong Zhang By placing code between these two statements, computations can be 13860f0b76eSHong Zhang done while messages are in transition. */ 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 14160f0b76eSHong Zhang 14260f0b76eSHong Zhang /* Get pointers to vector data */ 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localX,&x)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,mctx->U,&u)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Y,&y)); 14660f0b76eSHong Zhang 14760f0b76eSHong Zhang /* Get local grid boundaries */ 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 14960f0b76eSHong Zhang 15060f0b76eSHong Zhang /* Compute function over the locally owned part of the grid */ 15160f0b76eSHong Zhang for (j=ys; j<ys+ym; j++) { 15260f0b76eSHong Zhang for (i=xs; i<xs+xm; i++) { 15360f0b76eSHong Zhang uc = u[j][i].u; 15460f0b76eSHong Zhang ucb = x[j][i].u; 15560f0b76eSHong Zhang uxx = (-2.0*ucb + x[j][i-1].u + x[j][i+1].u)*sx; 15660f0b76eSHong Zhang uyy = (-2.0*ucb + x[j-1][i].u + x[j+1][i].u)*sy; 15760f0b76eSHong Zhang vc = u[j][i].v; 15860f0b76eSHong Zhang vcb = x[j][i].v; 15960f0b76eSHong Zhang vxx = (-2.0*vcb + x[j][i-1].v + x[j][i+1].v)*sx; 16060f0b76eSHong Zhang vyy = (-2.0*vcb + x[j-1][i].v + x[j+1][i].v)*sy; 16160f0b76eSHong Zhang y[j][i].u = appctx->D1*(uxx + uyy) - ucb*(vc*vc+appctx->gamma) + vc*vc*vcb; 16260f0b76eSHong Zhang y[j][i].v = appctx->D2*(vxx + vyy) - 2.0*uc*vc*ucb + (2.0*uc*vc-appctx->gamma - appctx->kappa)*vcb; 16360f0b76eSHong Zhang y[j][i].u = mctx->shift*ucb - y[j][i].u; 16460f0b76eSHong Zhang y[j][i].v = mctx->shift*vcb - y[j][i].v; 16560f0b76eSHong Zhang } 16660f0b76eSHong Zhang } 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,mctx->U,&u)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Y,&y)); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 17160f0b76eSHong Zhang PetscFunctionReturn(0); 17260f0b76eSHong Zhang } 17360f0b76eSHong Zhang 17460f0b76eSHong Zhang static PetscErrorCode MyIMatMult(Mat A_shell,Vec X,Vec Y) 17560f0b76eSHong Zhang { 17660f0b76eSHong Zhang MCtx *mctx; 17760f0b76eSHong Zhang AppCtx *appctx; 17860f0b76eSHong Zhang DM da; 17960f0b76eSHong Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 18060f0b76eSHong Zhang PetscReal hx,hy,sx,sy; 18160f0b76eSHong Zhang PetscScalar uc,uxx,uyy,vc,vxx,vyy,ucb,vcb; 18260f0b76eSHong Zhang Field **u,**x,**y; 18360f0b76eSHong Zhang Vec localX; 18460f0b76eSHong Zhang 18560f0b76eSHong Zhang PetscFunctionBeginUser; 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A_shell,&mctx)); 18760f0b76eSHong Zhang appctx = mctx->appctx; 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(mctx->ts,&da)); 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 190*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)); 19160f0b76eSHong Zhang hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 19260f0b76eSHong Zhang hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 19360f0b76eSHong Zhang 19460f0b76eSHong Zhang /* Scatter ghost points to local vector,using the 2-step process 19560f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 19660f0b76eSHong Zhang By placing code between these two statements, computations can be 19760f0b76eSHong Zhang done while messages are in transition. */ 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 20060f0b76eSHong Zhang 20160f0b76eSHong Zhang /* Get pointers to vector data */ 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localX,&x)); 203*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,mctx->U,&u)); 204*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,Y,&y)); 20560f0b76eSHong Zhang 20660f0b76eSHong Zhang /* Get local grid boundaries */ 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 20860f0b76eSHong Zhang 20960f0b76eSHong Zhang /* Compute function over the locally owned part of the grid */ 21060f0b76eSHong Zhang for (j=ys; j<ys+ym; j++) { 21160f0b76eSHong Zhang for (i=xs; i<xs+xm; i++) { 21260f0b76eSHong Zhang uc = u[j][i].u; 21360f0b76eSHong Zhang ucb = x[j][i].u; 21460f0b76eSHong Zhang uxx = (-2.0*ucb + x[j][i-1].u + x[j][i+1].u)*sx; 21560f0b76eSHong Zhang uyy = (-2.0*ucb + x[j-1][i].u + x[j+1][i].u)*sy; 21660f0b76eSHong Zhang vc = u[j][i].v; 21760f0b76eSHong Zhang vcb = x[j][i].v; 21860f0b76eSHong Zhang vxx = (-2.0*vcb + x[j][i-1].v + x[j][i+1].v)*sx; 21960f0b76eSHong Zhang vyy = (-2.0*vcb + x[j-1][i].v + x[j+1][i].v)*sy; 22060f0b76eSHong Zhang y[j][i].u = appctx->D1*(uxx + uyy) - (vc*vc+appctx->gamma)*ucb - 2.0*uc*vc*vcb; 22160f0b76eSHong Zhang y[j][i].v = appctx->D2*(vxx + vyy) + vc*vc*ucb + (2.0*uc*vc-appctx->gamma-appctx->kappa)*vcb; 22260f0b76eSHong Zhang y[j][i].u = mctx->shift*ucb - y[j][i].u; 22360f0b76eSHong Zhang y[j][i].v = mctx->shift*vcb - y[j][i].v; 22460f0b76eSHong Zhang } 22560f0b76eSHong Zhang } 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,mctx->U,&u)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,Y,&y)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 23060f0b76eSHong Zhang PetscFunctionReturn(0); 23160f0b76eSHong Zhang } 23260f0b76eSHong Zhang 23360f0b76eSHong Zhang int main(int argc,char **argv) 23460f0b76eSHong Zhang { 23560f0b76eSHong Zhang TS ts; /* ODE integrator */ 23660f0b76eSHong Zhang Vec x; /* solution */ 23760f0b76eSHong Zhang PetscErrorCode ierr; 23860f0b76eSHong Zhang DM da; 23960f0b76eSHong Zhang AppCtx appctx; 24060f0b76eSHong Zhang MCtx mctx; 24160f0b76eSHong Zhang Vec lambda[1]; 24260f0b76eSHong Zhang PetscBool forwardonly=PETSC_FALSE,implicitform=PETSC_TRUE,mf = PETSC_FALSE; 24360f0b76eSHong Zhang PetscLogDouble v1,v2; 24460f0b76eSHong Zhang #if defined(PETSC_USE_LOG) 24560f0b76eSHong Zhang PetscLogStage stage; 24660f0b76eSHong Zhang #endif 24760f0b76eSHong Zhang 24860f0b76eSHong Zhang ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL)); 25160f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE; 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL)); 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-mf",&mf,NULL)); 25460f0b76eSHong Zhang 25560f0b76eSHong Zhang appctx.D1 = 8.0e-5; 25660f0b76eSHong Zhang appctx.D2 = 4.0e-5; 25760f0b76eSHong Zhang appctx.gamma = .024; 25860f0b76eSHong Zhang appctx.kappa = .06; 25960f0b76eSHong Zhang 26060f0b76eSHong Zhang PetscLogStageRegister("MyAdjoint", &stage); 26160f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 26260f0b76eSHong Zhang Create distributed array (DMDA) to manage parallel grid and vectors 26360f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da)); 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 266*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 267*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"u")); 268*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,1,"v")); 26960f0b76eSHong Zhang 27060f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27160f0b76eSHong Zhang Extract global vectors from DMDA; then duplicate for remaining 27260f0b76eSHong Zhang vectors that are the same types 27360f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 27560f0b76eSHong Zhang 27660f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27760f0b76eSHong Zhang Create timestepping solver context 27860f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 279*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 281*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 28360f0b76eSHong Zhang if (!implicitform) { 284*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSRK)); 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&appctx)); 286*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx)); 28760f0b76eSHong Zhang } else { 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSCN)); 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&appctx)); 29060f0b76eSHong Zhang if (appctx.aijpc) { 29160f0b76eSHong Zhang Mat A,B; 29260f0b76eSHong Zhang 293*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATSELL)); 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&A)); 295*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 29660f0b76eSHong Zhang /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */ 297*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,A,B,IJacobian,&appctx)); 298*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 299*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 30060f0b76eSHong Zhang } else { 301*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,NULL,NULL,IJacobian,&appctx)); 30260f0b76eSHong Zhang } 30360f0b76eSHong Zhang } 30460f0b76eSHong Zhang 30560f0b76eSHong Zhang if (mf) { 30660f0b76eSHong Zhang PetscInt xm,ym,Mx,My,dof; 30760f0b76eSHong Zhang mctx.ts = ts; 30860f0b76eSHong Zhang mctx.appctx = &appctx; 309*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&mctx.U)); 31060f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 31160f0b76eSHong Zhang Create matrix free context 31260f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,&dof,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,NULL,NULL,NULL,&xm,&ym,NULL)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,dof*xm*ym,PETSC_DETERMINE,dof*Mx*My,dof*Mx*My,&mctx,&appctx.A)); 316*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(appctx.A,MATOP_MULT_TRANSPOSE,(void (*)(void))MyRHSMatMultTranspose)); 31760f0b76eSHong Zhang if (!implicitform) { /* for explicit methods only */ 318*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,appctx.A,appctx.A,RHSJacobianShell,&appctx)); 31960f0b76eSHong Zhang } else { 320*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecDuplicate(appctx.U,&mctx.Udot)); */ 321*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(appctx.A,MATOP_MULT,(void (*)(void))MyIMatMult)); 322*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(appctx.A,MATOP_MULT_TRANSPOSE,(void (*)(void))MyIMatMultTranspose)); 323*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,appctx.A,appctx.A,IJacobianShell,&appctx)); 32460f0b76eSHong Zhang } 32560f0b76eSHong Zhang } 32660f0b76eSHong Zhang 32760f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 32860f0b76eSHong Zhang Set initial conditions 32960f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(da,x)); 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,x)); 33260f0b76eSHong Zhang 33360f0b76eSHong Zhang /* 33460f0b76eSHong Zhang Have the TS save its trajectory so that TSAdjointSolve() may be used 33560f0b76eSHong Zhang */ 336*5f80ce2aSJacob Faibussowitsch if (!forwardonly) CHKERRQ(TSSetSaveTrajectory(ts)); 33760f0b76eSHong Zhang 33860f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 33960f0b76eSHong Zhang Set solver options 34060f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 341*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,200.0)); 342*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.5)); 343*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 344*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 34560f0b76eSHong Zhang 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTime(&v1)); 34760f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 34860f0b76eSHong Zhang Solve ODE system 34960f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,x)); 35160f0b76eSHong Zhang if (!forwardonly) { 35260f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 35360f0b76eSHong Zhang Start the Adjoint model 35460f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 355*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&lambda[0])); 35660f0b76eSHong Zhang /* Reset initial conditions for the adjoint integration */ 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(InitializeLambda(da,lambda[0],0.5,0.5)); 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(ts,1,lambda,NULL)); 35960f0b76eSHong Zhang PetscLogStagePush(stage); 360*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(ts)); 36160f0b76eSHong Zhang PetscLogStagePop(); 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&lambda[0])); 36360f0b76eSHong Zhang } 364*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTime(&v2)); 36560f0b76eSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD," %.3lf ",v2-v1); 36660f0b76eSHong Zhang 36760f0b76eSHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 36860f0b76eSHong Zhang Free work space. All PETSc objects should be destroyed when they 36960f0b76eSHong Zhang are no longer needed. 37060f0b76eSHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 371*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 372*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 373*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 37460f0b76eSHong Zhang if (mf) { 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&mctx.U)); 376*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecDestroy(&mctx.Udot));*/ 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&appctx.A)); 37860f0b76eSHong Zhang } 37960f0b76eSHong Zhang ierr = PetscFinalize(); 38060f0b76eSHong Zhang return ierr; 38160f0b76eSHong Zhang } 38260f0b76eSHong Zhang 38360f0b76eSHong Zhang /* ------------------------------------------------------------------- */ 38460f0b76eSHong Zhang PetscErrorCode RHSJacobianShell(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx) 38560f0b76eSHong Zhang { 38660f0b76eSHong Zhang MCtx *mctx; 38760f0b76eSHong Zhang 38860f0b76eSHong Zhang PetscFunctionBegin; 389*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&mctx)); 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,mctx->U)); 39160f0b76eSHong Zhang PetscFunctionReturn(0); 39260f0b76eSHong Zhang } 39360f0b76eSHong Zhang 39460f0b76eSHong Zhang PetscErrorCode IJacobianShell(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx) 39560f0b76eSHong Zhang { 39660f0b76eSHong Zhang MCtx *mctx; 39760f0b76eSHong Zhang 39860f0b76eSHong Zhang PetscFunctionBegin; 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(A,&mctx)); 400*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,mctx->U)); 401*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecCopy(Udot,mctx->Udot)); */ 40260f0b76eSHong Zhang mctx->shift = a; 40360f0b76eSHong Zhang PetscFunctionReturn(0); 40460f0b76eSHong Zhang } 40560f0b76eSHong Zhang 40660f0b76eSHong Zhang /* ------------------------------------------------------------------- */ 40760f0b76eSHong Zhang PetscErrorCode InitialConditions(DM da,Vec U) 40860f0b76eSHong Zhang { 40960f0b76eSHong Zhang PetscInt i,j,xs,ys,xm,ym,Mx,My; 41060f0b76eSHong Zhang Field **u; 41160f0b76eSHong Zhang PetscReal hx,hy,x,y; 41260f0b76eSHong Zhang 41360f0b76eSHong Zhang PetscFunctionBegin; 414*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)); 41560f0b76eSHong Zhang 41660f0b76eSHong Zhang hx = 2.5/(PetscReal)Mx; 41760f0b76eSHong Zhang hy = 2.5/(PetscReal)My; 41860f0b76eSHong Zhang 41960f0b76eSHong Zhang /* 42060f0b76eSHong Zhang Get pointers to vector data 42160f0b76eSHong Zhang */ 422*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 42360f0b76eSHong Zhang 42460f0b76eSHong Zhang /* 42560f0b76eSHong Zhang Get local grid boundaries 42660f0b76eSHong Zhang */ 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 42860f0b76eSHong Zhang 42960f0b76eSHong Zhang /* 43060f0b76eSHong Zhang Compute function over the locally owned part of the grid 43160f0b76eSHong Zhang */ 43260f0b76eSHong Zhang for (j=ys; j<ys+ym; j++) { 43360f0b76eSHong Zhang y = j*hy; 43460f0b76eSHong Zhang for (i=xs; i<xs+xm; i++) { 43560f0b76eSHong Zhang x = i*hx; 43660f0b76eSHong Zhang 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); 43760f0b76eSHong Zhang else u[j][i].v = 0.0; 43860f0b76eSHong Zhang 43960f0b76eSHong Zhang u[j][i].u = 1.0 - 2.0*u[j][i].v; 44060f0b76eSHong Zhang } 44160f0b76eSHong Zhang } 44260f0b76eSHong Zhang 44360f0b76eSHong Zhang /* 44460f0b76eSHong Zhang Restore vectors 44560f0b76eSHong Zhang */ 446*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 44760f0b76eSHong Zhang PetscFunctionReturn(0); 44860f0b76eSHong Zhang } 44960f0b76eSHong Zhang 45060f0b76eSHong Zhang /*TEST 45160f0b76eSHong Zhang 45260f0b76eSHong Zhang build: 45360f0b76eSHong Zhang depends: reaction_diffusion.c 45460f0b76eSHong Zhang requires: !complex !single 45560f0b76eSHong Zhang 45660f0b76eSHong Zhang test: 45760f0b76eSHong Zhang requires: cams 458533251d3SHong Zhang args: -mf -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -ts_trajectory_max_units_ram 6 -ts_trajectory_memory_type cams 45960f0b76eSHong Zhang TEST*/ 460