xref: /petsc/src/ts/tests/ex14.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(X,&n));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&xx));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&xx));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
315f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
535f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetType(ts,&type));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&U));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(U));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,Tf));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,NULL));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRollBack(ts));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,NULL));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&U));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(U,&X));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(TSEvaluateStep(ts,order,X,NULL));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&xx));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
76c4762a1bSJed Brown 
775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&U));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(U,&Y));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSEvaluateStep(ts,order-1,Y,&done));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Y,&yy));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Y));
87c4762a1bSJed Brown 
885f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts,&adapt));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdaptSetCheckStage(adapt,TestCheckStage));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetErrorIfStepFails(ts,PETSC_FALSE));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,NULL));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdaptSetCheckStage(adapt,NULL));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetErrorIfStepFails(ts,PETSC_TRUE));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetConvergedReason(ts,&reason));
993c633725SBarry Smith   PetscCheck(step == 2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bad step number %D, %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;
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRKSetType(ts,type));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRKGetType(ts,&type));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRKGetOrder(ts,&order));
116c4762a1bSJed Brown 
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts,&adapt));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdaptGetType(adapt,&adapttype));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrncpy(savetype,adapttype,sizeof(savetype)));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(type,TSRK1FE,&rk1));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(type,TSRK3,&rk3));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(type,TSRK4,&rk4));
1235f80ce2aSJacob Faibussowitsch   if (rk1 || rk3 || rk4) CHKERRQ(TSAdaptSetType(adapt,TSADAPTNONE));
124c4762a1bSJed Brown 
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(TestExplicitTS(ts,order,type));
126c4762a1bSJed Brown 
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts,&adapt));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
138*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
139c4762a1bSJed Brown 
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSRK));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N,&X));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,X));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
147c4762a1bSJed Brown 
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK1FE));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK2A));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK3));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK3BS));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK4));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK5F));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK5DP));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK5BS));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK6VR));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK7VR));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(TestTSRK(ts,TSRK8VR));
159c4762a1bSJed Brown 
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRollBack(ts));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
162c4762a1bSJed Brown 
163*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
164*b122ec5aSJacob 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