static char help[] = "Test TSComputeIJacobian()\n\n"; #include #include #include #include typedef struct { PetscScalar u,v; } Field; typedef struct { PetscReal D1,D2,gamma,kappa; } AppCtx; PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx) { AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ DM da; PetscInt i,j,Mx,My,xs,ys,xm,ym; PetscReal hx,hy,sx,sy; PetscScalar uc,vc; Field **u; Vec localU; MatStencil stencil[6],rowstencil; PetscScalar entries[6]; PetscFunctionBegin; PetscCall(TSGetDM(ts,&da)); PetscCall(DMGetLocalVector(da,&localU)); PetscCall(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)); hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); /* Scatter ghost points to local vector,using the 2-step process DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). By placing code between these two statements, computations can be done while messages are in transition. */ PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); /* Get pointers to vector data */ PetscCall(DMDAVecGetArrayRead(da,localU,&u)); /* Get local grid boundaries */ PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); stencil[0].k = 0; stencil[1].k = 0; stencil[2].k = 0; stencil[3].k = 0; stencil[4].k = 0; stencil[5].k = 0; rowstencil.k = 0; /* Compute function over the locally owned part of the grid */ for (j=ys; jD1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy; stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy; stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx; stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx; stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma; stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc; rowstencil.i = i; rowstencil.c = 0; PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES)); stencil[0].c = 1; entries[0] = appctx->D2*sy; stencil[1].c = 1; entries[1] = appctx->D2*sy; stencil[2].c = 1; entries[2] = appctx->D2*sx; stencil[3].c = 1; entries[3] = appctx->D2*sx; stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa; stencil[5].c = 0; entries[5] = vc*vc; rowstencil.c = 1; PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES)); /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ } } /* Restore vectors */ PetscCall(PetscLogFlops(19*xm*ym)); PetscCall(DMDAVecRestoreArrayRead(da,localU,&u)); PetscCall(DMRestoreLocalVector(da,&localU)); PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); PetscFunctionReturn(0); } PetscErrorCode InitialConditions(DM da,Vec U) { PetscInt i,j,xs,ys,xm,ym,Mx,My; Field **u; PetscReal hx,hy,x,y; PetscFunctionBegin; PetscCall(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)); hx = 2.5/(PetscReal)Mx; hy = 2.5/(PetscReal)My; /* Get pointers to vector data */ PetscCall(DMDAVecGetArray(da,U,&u)); /* Get local grid boundaries */ PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); /* Compute function over the locally owned part of the grid */ for (j=ys; j 100.0*PETSC_MACHINE_EPSILON) { PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error Norm %g \n Incorrect behaviour of TSComputeIJacobian(). The two matrices should have the same results.\n",norm)); } PetscCall(MatDestroy(&Jac)); PetscCall(MatDestroy(&Jac2)); PetscCall(VecDestroy(&U)); PetscCall(VecDestroy(&Udot)); PetscCall(TSDestroy(&ts)); PetscCall(DMDestroy(&da)); PetscCall(PetscFinalize()); return 0; } /*TEST test: TEST*/