15591c571SHong Zhang static char help[] = "Test TSComputeIJacobian()\n\n"; 25591c571SHong Zhang 35591c571SHong Zhang #include <petscsys.h> 45591c571SHong Zhang #include <petscdm.h> 55591c571SHong Zhang #include <petscdmda.h> 65591c571SHong Zhang #include <petscts.h> 75591c571SHong Zhang 85591c571SHong Zhang typedef struct { 95591c571SHong Zhang PetscScalar u,v; 105591c571SHong Zhang } Field; 115591c571SHong Zhang 125591c571SHong Zhang typedef struct { 135591c571SHong Zhang PetscReal D1,D2,gamma,kappa; 145591c571SHong Zhang } AppCtx; 155591c571SHong Zhang 165591c571SHong Zhang PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx) 175591c571SHong Zhang { 185591c571SHong Zhang AppCtx *appctx = (AppCtx*)ctx; /* user-defined application context */ 195591c571SHong Zhang DM da; 205591c571SHong Zhang PetscInt i,j,Mx,My,xs,ys,xm,ym; 215591c571SHong Zhang PetscReal hx,hy,sx,sy; 225591c571SHong Zhang PetscScalar uc,vc; 235591c571SHong Zhang Field **u; 245591c571SHong Zhang Vec localU; 255591c571SHong Zhang MatStencil stencil[6],rowstencil; 265591c571SHong Zhang PetscScalar entries[6]; 275591c571SHong Zhang 285591c571SHong Zhang PetscFunctionBegin; 295f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localU)); 315f80ce2aSJacob 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)); 325591c571SHong Zhang 335591c571SHong Zhang hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx); 345591c571SHong Zhang hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy); 355591c571SHong Zhang 365591c571SHong Zhang /* 375591c571SHong Zhang Scatter ghost points to local vector,using the 2-step process 385591c571SHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 395591c571SHong Zhang By placing code between these two statements, computations can be 405591c571SHong Zhang done while messages are in transition. 415591c571SHong Zhang */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 445591c571SHong Zhang 455591c571SHong Zhang /* 465591c571SHong Zhang Get pointers to vector data 475591c571SHong Zhang */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localU,&u)); 495591c571SHong Zhang 505591c571SHong Zhang /* 515591c571SHong Zhang Get local grid boundaries 525591c571SHong Zhang */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 545591c571SHong Zhang 555591c571SHong Zhang stencil[0].k = 0; 565591c571SHong Zhang stencil[1].k = 0; 575591c571SHong Zhang stencil[2].k = 0; 585591c571SHong Zhang stencil[3].k = 0; 595591c571SHong Zhang stencil[4].k = 0; 605591c571SHong Zhang stencil[5].k = 0; 615591c571SHong Zhang rowstencil.k = 0; 625591c571SHong Zhang /* 635591c571SHong Zhang Compute function over the locally owned part of the grid 645591c571SHong Zhang */ 655591c571SHong Zhang for (j=ys; j<ys+ym; j++) { 665591c571SHong Zhang 675591c571SHong Zhang stencil[0].j = j-1; 685591c571SHong Zhang stencil[1].j = j+1; 695591c571SHong Zhang stencil[2].j = j; 705591c571SHong Zhang stencil[3].j = j; 715591c571SHong Zhang stencil[4].j = j; 725591c571SHong Zhang stencil[5].j = j; 735591c571SHong Zhang rowstencil.k = 0; rowstencil.j = j; 745591c571SHong Zhang for (i=xs; i<xs+xm; i++) { 755591c571SHong Zhang uc = u[j][i].u; 765591c571SHong Zhang vc = u[j][i].v; 775591c571SHong Zhang 785591c571SHong Zhang /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 795591c571SHong Zhang uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 805591c571SHong Zhang 815591c571SHong Zhang vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 825591c571SHong Zhang vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 835591c571SHong Zhang f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 845591c571SHong Zhang 855591c571SHong Zhang stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy; 865591c571SHong Zhang stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy; 875591c571SHong Zhang stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx; 885591c571SHong Zhang stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx; 895591c571SHong Zhang stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma; 905591c571SHong Zhang stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc; 915591c571SHong Zhang rowstencil.i = i; rowstencil.c = 0; 925591c571SHong Zhang 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES)); 945591c571SHong Zhang stencil[0].c = 1; entries[0] = appctx->D2*sy; 955591c571SHong Zhang stencil[1].c = 1; entries[1] = appctx->D2*sy; 965591c571SHong Zhang stencil[2].c = 1; entries[2] = appctx->D2*sx; 975591c571SHong Zhang stencil[3].c = 1; entries[3] = appctx->D2*sx; 985591c571SHong Zhang stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa; 995591c571SHong Zhang stencil[5].c = 0; entries[5] = vc*vc; 1005591c571SHong Zhang rowstencil.c = 1; 1015591c571SHong Zhang 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES)); 1035591c571SHong Zhang /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 1045591c571SHong Zhang } 1055591c571SHong Zhang } 1065591c571SHong Zhang 1075591c571SHong Zhang /* 1085591c571SHong Zhang Restore vectors 1095591c571SHong Zhang */ 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogFlops(19*xm*ym)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localU)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 1165591c571SHong Zhang PetscFunctionReturn(0); 1175591c571SHong Zhang } 1185591c571SHong Zhang 1195591c571SHong Zhang PetscErrorCode InitialConditions(DM da,Vec U) 1205591c571SHong Zhang { 1215591c571SHong Zhang PetscInt i,j,xs,ys,xm,ym,Mx,My; 1225591c571SHong Zhang Field **u; 1235591c571SHong Zhang PetscReal hx,hy,x,y; 1245591c571SHong Zhang 1255591c571SHong Zhang PetscFunctionBegin; 1265f80ce2aSJacob 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)); 1275591c571SHong Zhang 1285591c571SHong Zhang hx = 2.5/(PetscReal)Mx; 1295591c571SHong Zhang hy = 2.5/(PetscReal)My; 1305591c571SHong Zhang 1315591c571SHong Zhang /* 1325591c571SHong Zhang Get pointers to vector data 1335591c571SHong Zhang */ 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 1355591c571SHong Zhang 1365591c571SHong Zhang /* 1375591c571SHong Zhang Get local grid boundaries 1385591c571SHong Zhang */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL)); 1405591c571SHong Zhang 1415591c571SHong Zhang /* 1425591c571SHong Zhang Compute function over the locally owned part of the grid 1435591c571SHong Zhang */ 1445591c571SHong Zhang for (j=ys; j<ys+ym; j++) { 1455591c571SHong Zhang y = j*hy; 1465591c571SHong Zhang for (i=xs; i<xs+xm; i++) { 1475591c571SHong Zhang x = i*hx; 1485591c571SHong 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); 1495591c571SHong Zhang else u[j][i].v = 0.0; 1505591c571SHong Zhang 1515591c571SHong Zhang u[j][i].u = 1.0 - 2.0*u[j][i].v; 1525591c571SHong Zhang } 1535591c571SHong Zhang } 1545591c571SHong Zhang 1555591c571SHong Zhang /* 1565591c571SHong Zhang Restore vectors 1575591c571SHong Zhang */ 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 1595591c571SHong Zhang PetscFunctionReturn(0); 1605591c571SHong Zhang } 1615591c571SHong Zhang 1625591c571SHong Zhang int main(int argc,char **argv) 1635591c571SHong Zhang { 1645591c571SHong Zhang TS ts; 1655591c571SHong Zhang Vec U,Udot; 1665591c571SHong Zhang Mat Jac,Jac2; 1675591c571SHong Zhang DM da; 1685591c571SHong Zhang AppCtx appctx; 1695591c571SHong Zhang PetscReal t = 0,shift,norm; 1705591c571SHong Zhang 171*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 1725591c571SHong Zhang 1735591c571SHong Zhang appctx.D1 = 8.0e-5; 1745591c571SHong Zhang appctx.D2 = 4.0e-5; 1755591c571SHong Zhang appctx.gamma = .024; 1765591c571SHong Zhang appctx.kappa = .06; 1775591c571SHong Zhang 1785591c571SHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1795591c571SHong Zhang Create distributed array (DMDA) to manage parallel grid and vectors 1805591c571SHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1815f80ce2aSJacob 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)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"u")); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,1,"v")); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(da,MATAIJ)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&Jac)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(da,&Jac2)); 1895591c571SHong Zhang 1905591c571SHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1915591c571SHong Zhang Extract global vectors from DMDA; then duplicate for remaining 1925591c571SHong Zhang vectors that are the same types 1935591c571SHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&U)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Udot)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(Udot,0.0)); 1975591c571SHong Zhang 1985591c571SHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1995591c571SHong Zhang Create timestepping solver context 2005591c571SHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,Jac,Jac,RHSJacobian,&appctx)); 2055591c571SHong Zhang 2065591c571SHong Zhang /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 2075591c571SHong Zhang Set initial conditions 2085591c571SHong Zhang - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(InitialConditions(da,U)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,U)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetUp(ts)); 2135591c571SHong Zhang 2145591c571SHong Zhang shift = 2.; 2155f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeIJacobian(ts,t,U,Udot,shift,Jac2,Jac2,PETSC_FALSE)); 2165591c571SHong Zhang shift = 1.; 2175f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeIJacobian(ts,t,U,Udot,shift,Jac,Jac,PETSC_FALSE)); 2185591c571SHong Zhang shift = 2.; 2195f80ce2aSJacob Faibussowitsch CHKERRQ(TSComputeIJacobian(ts,t,U,Udot,shift,Jac,Jac,PETSC_FALSE)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Jac,-1,Jac2,SAME_NONZERO_PATTERN)); 2215f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Jac,NORM_INFINITY,&norm)); 2225591c571SHong Zhang if (norm > 100.0*PETSC_MACHINE_EPSILON) { 2235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error Norm %g \n Incorrect behaviour of TSComputeIJacobian(). The two matrices should have the same results.\n",norm)); 2245591c571SHong Zhang } 2255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Jac)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Jac2)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Udot)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 231*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 232*b122ec5aSJacob Faibussowitsch return 0; 2335591c571SHong Zhang } 2345591c571SHong Zhang 2355591c571SHong Zhang /*TEST 2365591c571SHong Zhang test: 2375591c571SHong Zhang 2385591c571SHong Zhang TEST*/ 239