xref: /petsc/src/ts/tests/ex14.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 
11c4762a1bSJed Brown   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X,&n));
139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&xx));
149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&ff));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   if (n >= 1)
17c4762a1bSJed Brown     ff[0] = 1;
18c4762a1bSJed Brown   for (i = 1; i < n; i++)
19c4762a1bSJed Brown     ff[i] = (i+1)*(xx[i-1]+PetscPowReal(t,i))/2;
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&xx));
229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&ff));
23c4762a1bSJed Brown   PetscFunctionReturn(0);
24c4762a1bSJed Brown }
25c4762a1bSJed Brown 
26c4762a1bSJed Brown PetscErrorCode TestCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec X,PetscBool *accept)
27c4762a1bSJed Brown {
28c4762a1bSJed Brown   PetscInt       step;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
32c4762a1bSJed Brown   *accept = (step >= 2) ? PETSC_FALSE : PETSC_TRUE;
33c4762a1bSJed Brown   PetscFunctionReturn(0);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown static PetscErrorCode TestExplicitTS(TS ts,PetscInt order,const char subtype[])
37c4762a1bSJed Brown {
38c4762a1bSJed Brown   PetscInt          i;
39c4762a1bSJed Brown   PetscReal         t;
40c4762a1bSJed Brown   Vec               U,X,Y;
41c4762a1bSJed Brown   TSType            type;
42c4762a1bSJed Brown   PetscBool         done;
43c4762a1bSJed Brown   const PetscScalar *xx;
44c4762a1bSJed Brown   const PetscScalar *yy;
45c4762a1bSJed Brown   const PetscReal   Tf  = 1;
46c4762a1bSJed Brown   const PetscReal   dt  = Tf/8;
47c4762a1bSJed Brown   const PetscReal   eps = 100*PETSC_MACHINE_EPSILON;
48c4762a1bSJed Brown   TSAdapt           adapt;
49c4762a1bSJed Brown   PetscInt          step;
50c4762a1bSJed Brown   TSConvergedReason reason;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(TSGetType(ts,&type));
549566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
559566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(U));
569566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
579566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0));
589566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
599566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,Tf));
609566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
619566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,NULL));
629566063dSJacob Faibussowitsch   PetscCall(TSRollBack(ts));
639566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,NULL));
649566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&X));
689566063dSJacob Faibussowitsch   PetscCall(TSEvaluateStep(ts,order,X,NULL));
699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&xx));
70c4762a1bSJed Brown   for (i = 0; i < order; i++) {
71c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(xx[i]) - PetscPowReal(t,i+1));
723c633725SBarry Smith     PetscCheck(error <= eps,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad solution, error %g, %s '%s'",(double)error,type,subtype);
73c4762a1bSJed Brown   }
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&xx));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
789566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&Y));
799566063dSJacob Faibussowitsch   PetscCall(TSEvaluateStep(ts,order-1,Y,&done));
809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y,&yy));
81c4762a1bSJed Brown   for (i = 0; done && i < order-1; i++) {
82c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(yy[i]) - PetscPowReal(t,i+1));
833c633725SBarry Smith     PetscCheck(error <= eps,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad estimator, error %g, %s '%s'",(double)error,type,subtype);
84c4762a1bSJed Brown   }
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y,&yy));
869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
899566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetCheckStage(adapt,TestCheckStage));
909566063dSJacob Faibussowitsch   PetscCall(TSSetErrorIfStepFails(ts,PETSC_FALSE));
919566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
929566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0));
939566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
949566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,NULL));
959566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetCheckStage(adapt,NULL));
969566063dSJacob Faibussowitsch   PetscCall(TSSetErrorIfStepFails(ts,PETSC_TRUE));
979566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
989566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts,&reason));
99*63a3b9bcSJacob Faibussowitsch   PetscCheck(step == 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad step number %" PetscInt_FMT ", %s '%s'",step,type,subtype);
1003c633725SBarry Smith   PetscCheck(reason == TS_DIVERGED_STEP_REJECTED,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad reason %s, %s '%s'",TSConvergedReasons[reason],type,subtype);
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown static PetscErrorCode TestTSRK(TS ts,TSRKType type)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   PetscInt       order;
107c4762a1bSJed Brown   TSAdapt        adapt;
108c4762a1bSJed Brown   PetscBool      rk1,rk3,rk4;
109c4762a1bSJed Brown   TSAdaptType    adapttype;
110c4762a1bSJed Brown   char           savetype[32];
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   PetscFunctionBegin;
1139566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts,type));
1149566063dSJacob Faibussowitsch   PetscCall(TSRKGetType(ts,&type));
1159566063dSJacob Faibussowitsch   PetscCall(TSRKGetOrder(ts,&order));
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
1189566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetType(adapt,&adapttype));
1199566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(savetype,adapttype,sizeof(savetype)));
1209566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type,TSRK1FE,&rk1));
1219566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type,TSRK3,&rk3));
1229566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type,TSRK4,&rk4));
1239566063dSJacob Faibussowitsch   if (rk1 || rk3 || rk4) PetscCall(TSAdaptSetType(adapt,TSADAPTNONE));
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(TestExplicitTS(ts,order,type));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
1289566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt,savetype));
129c4762a1bSJed Brown   PetscFunctionReturn(0);
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown int main(int argc, char *argv[])
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   TS             ts;
135c4762a1bSJed Brown   Vec            X;
136c4762a1bSJed Brown   PetscInt       N = 9;
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
1419566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSRK));
1429566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
1439566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&X));
1449566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
1459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1469566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK1FE));
1499566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK2A));
1509566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK3));
1519566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK3BS));
1529566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK4));
1539566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK5F));
1549566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK5DP));
1559566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK5BS));
1569566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK6VR));
1579566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK7VR));
1589566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts,TSRK8VR));
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(TSRollBack(ts));
1619566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
162c4762a1bSJed Brown 
1639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
164b122ec5aSJacob Faibussowitsch   return 0;
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown /*TEST
168c4762a1bSJed Brown 
169c4762a1bSJed Brown testset:
170c4762a1bSJed Brown   output_file: output/ex14.out
171c4762a1bSJed Brown   test:
172c4762a1bSJed Brown     suffix: 0
173c4762a1bSJed Brown   test:
174c4762a1bSJed Brown     suffix: 1
175c4762a1bSJed Brown     args: -ts_adapt_type none
176c4762a1bSJed Brown   test:
177c4762a1bSJed Brown     suffix: 2
178c4762a1bSJed Brown     args: -ts_adapt_type basic
179c4762a1bSJed Brown   test:
180c4762a1bSJed Brown     suffix: 3
181c4762a1bSJed Brown     args: -ts_adapt_type dsp
182c4762a1bSJed Brown 
183c4762a1bSJed Brown TEST*/
184