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