xref: /petsc/src/ts/tests/ex14.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] ="Tests all TSRK types \n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscts.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   PetscInt          i,n;
8*c4762a1bSJed Brown   const PetscScalar *xx;
9*c4762a1bSJed Brown   /* */ PetscScalar *ff;
10*c4762a1bSJed Brown   PetscErrorCode    ierr;
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown   PetscFunctionBegin;
13*c4762a1bSJed Brown   ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr);
14*c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&xx);CHKERRQ(ierr);
15*c4762a1bSJed Brown   ierr = VecGetArray(F,&ff);CHKERRQ(ierr);
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown   if (n >= 1)
18*c4762a1bSJed Brown     ff[0] = 1;
19*c4762a1bSJed Brown   for (i = 1; i < n; i++)
20*c4762a1bSJed Brown     ff[i] = (i+1)*(xx[i-1]+PetscPowReal(t,i))/2;
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&xx);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&ff);CHKERRQ(ierr);
24*c4762a1bSJed Brown   PetscFunctionReturn(0);
25*c4762a1bSJed Brown }
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown PetscErrorCode TestCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec X,PetscBool *accept)
28*c4762a1bSJed Brown {
29*c4762a1bSJed Brown   PetscInt       step;
30*c4762a1bSJed Brown   PetscErrorCode ierr;
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown   PetscFunctionBegin;
33*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
34*c4762a1bSJed Brown   *accept = (step >= 2) ? PETSC_FALSE : PETSC_TRUE;
35*c4762a1bSJed Brown   PetscFunctionReturn(0);
36*c4762a1bSJed Brown }
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown static PetscErrorCode TestExplicitTS(TS ts,PetscInt order,const char subtype[])
39*c4762a1bSJed Brown {
40*c4762a1bSJed Brown   PetscInt          i;
41*c4762a1bSJed Brown   PetscReal         t;
42*c4762a1bSJed Brown   Vec               U,X,Y;
43*c4762a1bSJed Brown   TSType            type;
44*c4762a1bSJed Brown   PetscBool         done;
45*c4762a1bSJed Brown   const PetscScalar *xx;
46*c4762a1bSJed Brown   const PetscScalar *yy;
47*c4762a1bSJed Brown   const PetscReal   Tf  = 1;
48*c4762a1bSJed Brown   const PetscReal   dt  = Tf/8;
49*c4762a1bSJed Brown   const PetscReal   eps = 100*PETSC_MACHINE_EPSILON;
50*c4762a1bSJed Brown   TSAdapt           adapt;
51*c4762a1bSJed Brown   PetscInt          step;
52*c4762a1bSJed Brown   TSConvergedReason reason;
53*c4762a1bSJed Brown   PetscErrorCode    ierr;
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   PetscFunctionBegin;
56*c4762a1bSJed Brown   ierr = TSGetType(ts,&type);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = VecZeroEntries(U);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = TSSetTime(ts,0);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,Tf);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = TSSolve(ts,NULL);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = VecDuplicate(U,&X);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = TSEvaluateStep(ts,order,X,NULL);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&xx);CHKERRQ(ierr);
71*c4762a1bSJed Brown   for (i = 0; i < order; i++) {
72*c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(xx[i]) - PetscPowReal(t,i+1));
73*c4762a1bSJed Brown     if (error > eps) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad solution, error %g, %s '%s'",(double)error,type,subtype);
74*c4762a1bSJed Brown   }
75*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&xx);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = VecDuplicate(U,&Y);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = TSEvaluateStep(ts,order-1,Y,&done);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&yy);CHKERRQ(ierr);
82*c4762a1bSJed Brown   for (i = 0; done && i < order-1; i++) {
83*c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(yy[i]) - PetscPowReal(t,i+1));
84*c4762a1bSJed Brown     if (error > eps) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad estimator, error %g, %s '%s'",(double)error,type,subtype);
85*c4762a1bSJed Brown   }
86*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&yy);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = VecDestroy(&Y);CHKERRQ(ierr);
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = TSAdaptSetCheckStage(adapt,TestCheckStage);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = TSSetErrorIfStepFails(ts,PETSC_FALSE);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = TSSetTime(ts,0);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = TSSolve(ts,NULL);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = TSAdaptSetCheckStage(adapt,NULL);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = TSSetErrorIfStepFails(ts,PETSC_TRUE);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
100*c4762a1bSJed Brown   if (step != 2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad step number %D, %s '%s'",step,type,subtype);
101*c4762a1bSJed Brown   if (reason != TS_DIVERGED_STEP_REJECTED) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad reason %s, %s '%s'",TSConvergedReasons[reason],type,subtype);
102*c4762a1bSJed Brown   PetscFunctionReturn(0);
103*c4762a1bSJed Brown }
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown static PetscErrorCode TestTSRK(TS ts,TSRKType type)
106*c4762a1bSJed Brown {
107*c4762a1bSJed Brown   PetscInt       order;
108*c4762a1bSJed Brown   TSAdapt        adapt;
109*c4762a1bSJed Brown   PetscBool      rk1,rk3,rk4;
110*c4762a1bSJed Brown   TSAdaptType    adapttype;
111*c4762a1bSJed Brown   char           savetype[32];
112*c4762a1bSJed Brown   PetscErrorCode ierr;
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   PetscFunctionBegin;
115*c4762a1bSJed Brown   ierr = TSRKSetType(ts,type);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = TSRKGetType(ts,&type);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = TSRKGetOrder(ts,&order);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = TSAdaptGetType(adapt,&adapttype);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = PetscStrncpy(savetype,adapttype,sizeof(savetype));CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscStrcmp(type,TSRK1FE,&rk1);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = PetscStrcmp(type,TSRK3,&rk3);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = PetscStrcmp(type,TSRK4,&rk4);CHKERRQ(ierr);
125*c4762a1bSJed Brown   if (rk1 || rk3 || rk4) {ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);}
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   ierr = TestExplicitTS(ts,order,type);CHKERRQ(ierr);
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = TSAdaptSetType(adapt,savetype);CHKERRQ(ierr);
131*c4762a1bSJed Brown   PetscFunctionReturn(0);
132*c4762a1bSJed Brown }
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown int main(int argc, char *argv[])
135*c4762a1bSJed Brown {
136*c4762a1bSJed Brown   TS             ts;
137*c4762a1bSJed Brown   Vec            X;
138*c4762a1bSJed Brown   PetscInt       N = 9;
139*c4762a1bSJed Brown   PetscErrorCode ierr;
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,NULL);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,N,&X);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = VecDestroy(&X);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK1FE);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK2A);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK3);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK3BS);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK4);CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK5F);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK5DP);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK5BS);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK6VR);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK7VR);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = TestTSRK(ts,TSRK8VR);CHKERRQ(ierr);
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown   ierr = TSRollBack(ts);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   ierr = PetscFinalize();
167*c4762a1bSJed Brown   return ierr;
168*c4762a1bSJed Brown }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown /*TEST
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown testset:
173*c4762a1bSJed Brown   output_file: output/ex14.out
174*c4762a1bSJed Brown   test:
175*c4762a1bSJed Brown     suffix: 0
176*c4762a1bSJed Brown   test:
177*c4762a1bSJed Brown     suffix: 1
178*c4762a1bSJed Brown     args: -ts_adapt_type none
179*c4762a1bSJed Brown   test:
180*c4762a1bSJed Brown     suffix: 2
181*c4762a1bSJed Brown     args: -ts_adapt_type basic
182*c4762a1bSJed Brown   test:
183*c4762a1bSJed Brown     suffix: 3
184*c4762a1bSJed Brown     args: -ts_adapt_type dsp
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown TEST*/
187