xref: /petsc/src/ts/tests/ex24.c (revision 5591c571caaed1c8a2e812535ab4a129abe5ad43)
1*5591c571SHong Zhang static char help[] = "Test TSComputeIJacobian()\n\n";
2*5591c571SHong Zhang 
3*5591c571SHong Zhang #include <petscsys.h>
4*5591c571SHong Zhang #include <petscdm.h>
5*5591c571SHong Zhang #include <petscdmda.h>
6*5591c571SHong Zhang #include <petscts.h>
7*5591c571SHong Zhang 
8*5591c571SHong Zhang typedef struct {
9*5591c571SHong Zhang   PetscScalar u,v;
10*5591c571SHong Zhang } Field;
11*5591c571SHong Zhang 
12*5591c571SHong Zhang typedef struct {
13*5591c571SHong Zhang   PetscReal D1,D2,gamma,kappa;
14*5591c571SHong Zhang } AppCtx;
15*5591c571SHong Zhang 
16*5591c571SHong Zhang 
17*5591c571SHong Zhang PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx)
18*5591c571SHong Zhang {
19*5591c571SHong Zhang   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
20*5591c571SHong Zhang   DM             da;
21*5591c571SHong Zhang   PetscErrorCode ierr;
22*5591c571SHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
23*5591c571SHong Zhang   PetscReal      hx,hy,sx,sy;
24*5591c571SHong Zhang   PetscScalar    uc,vc;
25*5591c571SHong Zhang   Field          **u;
26*5591c571SHong Zhang   Vec            localU;
27*5591c571SHong Zhang   MatStencil     stencil[6],rowstencil;
28*5591c571SHong Zhang   PetscScalar    entries[6];
29*5591c571SHong Zhang 
30*5591c571SHong Zhang   PetscFunctionBegin;
31*5591c571SHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
32*5591c571SHong Zhang   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
33*5591c571SHong Zhang   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);
34*5591c571SHong Zhang 
35*5591c571SHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
36*5591c571SHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
37*5591c571SHong Zhang 
38*5591c571SHong Zhang   /*
39*5591c571SHong Zhang      Scatter ghost points to local vector,using the 2-step process
40*5591c571SHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
41*5591c571SHong Zhang      By placing code between these two statements, computations can be
42*5591c571SHong Zhang      done while messages are in transition.
43*5591c571SHong Zhang   */
44*5591c571SHong Zhang   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
45*5591c571SHong Zhang   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
46*5591c571SHong Zhang 
47*5591c571SHong Zhang   /*
48*5591c571SHong Zhang      Get pointers to vector data
49*5591c571SHong Zhang   */
50*5591c571SHong Zhang   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
51*5591c571SHong Zhang 
52*5591c571SHong Zhang   /*
53*5591c571SHong Zhang      Get local grid boundaries
54*5591c571SHong Zhang   */
55*5591c571SHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
56*5591c571SHong Zhang 
57*5591c571SHong Zhang   stencil[0].k = 0;
58*5591c571SHong Zhang   stencil[1].k = 0;
59*5591c571SHong Zhang   stencil[2].k = 0;
60*5591c571SHong Zhang   stencil[3].k = 0;
61*5591c571SHong Zhang   stencil[4].k = 0;
62*5591c571SHong Zhang   stencil[5].k = 0;
63*5591c571SHong Zhang   rowstencil.k = 0;
64*5591c571SHong Zhang   /*
65*5591c571SHong Zhang      Compute function over the locally owned part of the grid
66*5591c571SHong Zhang   */
67*5591c571SHong Zhang   for (j=ys; j<ys+ym; j++) {
68*5591c571SHong Zhang 
69*5591c571SHong Zhang     stencil[0].j = j-1;
70*5591c571SHong Zhang     stencil[1].j = j+1;
71*5591c571SHong Zhang     stencil[2].j = j;
72*5591c571SHong Zhang     stencil[3].j = j;
73*5591c571SHong Zhang     stencil[4].j = j;
74*5591c571SHong Zhang     stencil[5].j = j;
75*5591c571SHong Zhang     rowstencil.k = 0; rowstencil.j = j;
76*5591c571SHong Zhang     for (i=xs; i<xs+xm; i++) {
77*5591c571SHong Zhang       uc = u[j][i].u;
78*5591c571SHong Zhang       vc = u[j][i].v;
79*5591c571SHong Zhang 
80*5591c571SHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
81*5591c571SHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
82*5591c571SHong Zhang 
83*5591c571SHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
84*5591c571SHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
85*5591c571SHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
86*5591c571SHong Zhang 
87*5591c571SHong Zhang       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
88*5591c571SHong Zhang       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
89*5591c571SHong Zhang       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
90*5591c571SHong Zhang       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
91*5591c571SHong Zhang       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
92*5591c571SHong Zhang       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
93*5591c571SHong Zhang       rowstencil.i = i; rowstencil.c = 0;
94*5591c571SHong Zhang 
95*5591c571SHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
96*5591c571SHong Zhang       stencil[0].c = 1; entries[0] = appctx->D2*sy;
97*5591c571SHong Zhang       stencil[1].c = 1; entries[1] = appctx->D2*sy;
98*5591c571SHong Zhang       stencil[2].c = 1; entries[2] = appctx->D2*sx;
99*5591c571SHong Zhang       stencil[3].c = 1; entries[3] = appctx->D2*sx;
100*5591c571SHong Zhang       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
101*5591c571SHong Zhang       stencil[5].c = 0; entries[5] = vc*vc;
102*5591c571SHong Zhang       rowstencil.c = 1;
103*5591c571SHong Zhang 
104*5591c571SHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
105*5591c571SHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
106*5591c571SHong Zhang     }
107*5591c571SHong Zhang   }
108*5591c571SHong Zhang 
109*5591c571SHong Zhang   /*
110*5591c571SHong Zhang      Restore vectors
111*5591c571SHong Zhang   */
112*5591c571SHong Zhang   ierr = PetscLogFlops(19*xm*ym);CHKERRQ(ierr);
113*5591c571SHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
114*5591c571SHong Zhang   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
115*5591c571SHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116*5591c571SHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
117*5591c571SHong Zhang   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
118*5591c571SHong Zhang   PetscFunctionReturn(0);
119*5591c571SHong Zhang }
120*5591c571SHong Zhang 
121*5591c571SHong Zhang PetscErrorCode InitialConditions(DM da,Vec U)
122*5591c571SHong Zhang {
123*5591c571SHong Zhang   PetscErrorCode ierr;
124*5591c571SHong Zhang   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
125*5591c571SHong Zhang   Field          **u;
126*5591c571SHong Zhang   PetscReal      hx,hy,x,y;
127*5591c571SHong Zhang 
128*5591c571SHong Zhang   PetscFunctionBegin;
129*5591c571SHong Zhang   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);
130*5591c571SHong Zhang 
131*5591c571SHong Zhang   hx = 2.5/(PetscReal)Mx;
132*5591c571SHong Zhang   hy = 2.5/(PetscReal)My;
133*5591c571SHong Zhang 
134*5591c571SHong Zhang   /*
135*5591c571SHong Zhang      Get pointers to vector data
136*5591c571SHong Zhang   */
137*5591c571SHong Zhang   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
138*5591c571SHong Zhang 
139*5591c571SHong Zhang   /*
140*5591c571SHong Zhang      Get local grid boundaries
141*5591c571SHong Zhang   */
142*5591c571SHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
143*5591c571SHong Zhang 
144*5591c571SHong Zhang   /*
145*5591c571SHong Zhang      Compute function over the locally owned part of the grid
146*5591c571SHong Zhang   */
147*5591c571SHong Zhang   for (j=ys; j<ys+ym; j++) {
148*5591c571SHong Zhang     y = j*hy;
149*5591c571SHong Zhang     for (i=xs; i<xs+xm; i++) {
150*5591c571SHong Zhang       x = i*hx;
151*5591c571SHong 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);
152*5591c571SHong Zhang       else u[j][i].v = 0.0;
153*5591c571SHong Zhang 
154*5591c571SHong Zhang       u[j][i].u = 1.0 - 2.0*u[j][i].v;
155*5591c571SHong Zhang     }
156*5591c571SHong Zhang   }
157*5591c571SHong Zhang 
158*5591c571SHong Zhang   /*
159*5591c571SHong Zhang      Restore vectors
160*5591c571SHong Zhang   */
161*5591c571SHong Zhang   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
162*5591c571SHong Zhang   PetscFunctionReturn(0);
163*5591c571SHong Zhang }
164*5591c571SHong Zhang 
165*5591c571SHong Zhang int main(int argc,char **argv)
166*5591c571SHong Zhang {
167*5591c571SHong Zhang   TS             ts;
168*5591c571SHong Zhang   Vec            U,Udot;
169*5591c571SHong Zhang   Mat            Jac,Jac2;
170*5591c571SHong Zhang   DM             da;
171*5591c571SHong Zhang   AppCtx         appctx;
172*5591c571SHong Zhang   PetscReal      t = 0,shift,norm;
173*5591c571SHong Zhang   PetscErrorCode ierr;
174*5591c571SHong Zhang 
175*5591c571SHong Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
176*5591c571SHong Zhang 
177*5591c571SHong Zhang   appctx.D1    = 8.0e-5;
178*5591c571SHong Zhang   appctx.D2    = 4.0e-5;
179*5591c571SHong Zhang   appctx.gamma = .024;
180*5591c571SHong Zhang   appctx.kappa = .06;
181*5591c571SHong Zhang 
182*5591c571SHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183*5591c571SHong Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
184*5591c571SHong Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185*5591c571SHong Zhang   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);
186*5591c571SHong Zhang   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
187*5591c571SHong Zhang   ierr = DMSetUp(da);CHKERRQ(ierr);
188*5591c571SHong Zhang   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
189*5591c571SHong Zhang   ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr);
190*5591c571SHong Zhang   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
191*5591c571SHong Zhang   ierr = DMCreateMatrix(da,&Jac);CHKERRQ(ierr);
192*5591c571SHong Zhang   ierr = DMCreateMatrix(da,&Jac2);CHKERRQ(ierr);
193*5591c571SHong Zhang 
194*5591c571SHong Zhang   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195*5591c571SHong Zhang      Extract global vectors from DMDA; then duplicate for remaining
196*5591c571SHong Zhang      vectors that are the same types
197*5591c571SHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198*5591c571SHong Zhang   ierr = DMCreateGlobalVector(da,&U);CHKERRQ(ierr);
199*5591c571SHong Zhang   ierr = VecDuplicate(U,&Udot);CHKERRQ(ierr);
200*5591c571SHong Zhang   ierr = VecSet(Udot,0.0);CHKERRQ(ierr);
201*5591c571SHong Zhang 
202*5591c571SHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203*5591c571SHong Zhang      Create timestepping solver context
204*5591c571SHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205*5591c571SHong Zhang   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
206*5591c571SHong Zhang   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
207*5591c571SHong Zhang   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
208*5591c571SHong Zhang   ierr = TSSetRHSJacobian(ts,Jac,Jac,RHSJacobian,&appctx);CHKERRQ(ierr);
209*5591c571SHong Zhang 
210*5591c571SHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211*5591c571SHong Zhang      Set initial conditions
212*5591c571SHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213*5591c571SHong Zhang   ierr = InitialConditions(da,U);CHKERRQ(ierr);
214*5591c571SHong Zhang   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
215*5591c571SHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
216*5591c571SHong Zhang   ierr = TSSetUp(ts);CHKERRQ(ierr);
217*5591c571SHong Zhang 
218*5591c571SHong Zhang   shift = 2.;
219*5591c571SHong Zhang   ierr  = TSComputeIJacobian(ts,t,U,Udot,shift,Jac2,Jac2,PETSC_FALSE);CHKERRQ(ierr);
220*5591c571SHong Zhang   shift = 1.;
221*5591c571SHong Zhang   ierr  = TSComputeIJacobian(ts,t,U,Udot,shift,Jac,Jac,PETSC_FALSE);CHKERRQ(ierr);
222*5591c571SHong Zhang   shift = 2.;
223*5591c571SHong Zhang   ierr = TSComputeIJacobian(ts,t,U,Udot,shift,Jac,Jac,PETSC_FALSE);CHKERRQ(ierr);
224*5591c571SHong Zhang   ierr = MatAXPY(Jac,-1,Jac2,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
225*5591c571SHong Zhang   ierr = MatNorm(Jac,NORM_INFINITY,&norm);CHKERRQ(ierr);
226*5591c571SHong Zhang   if (norm > 100.0*PETSC_MACHINE_EPSILON) {
227*5591c571SHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error Norm %g \n Incorrect behaviour of TSComputeIJacobian(). The two matrices should have the same results.\n",norm);CHKERRQ(ierr);
228*5591c571SHong Zhang   }
229*5591c571SHong Zhang   ierr = MatDestroy(&Jac);CHKERRQ(ierr);
230*5591c571SHong Zhang   ierr = MatDestroy(&Jac2);CHKERRQ(ierr);
231*5591c571SHong Zhang   ierr = VecDestroy(&U);CHKERRQ(ierr);
232*5591c571SHong Zhang   ierr = VecDestroy(&Udot);CHKERRQ(ierr);
233*5591c571SHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
234*5591c571SHong Zhang   ierr = DMDestroy(&da);CHKERRQ(ierr);
235*5591c571SHong Zhang   ierr = PetscFinalize();
236*5591c571SHong Zhang   return(ierr);
237*5591c571SHong Zhang }
238*5591c571SHong Zhang 
239*5591c571SHong Zhang /*TEST
240*5591c571SHong Zhang   test:
241*5591c571SHong Zhang 
242*5591c571SHong Zhang TEST*/
243