xref: /petsc/src/ts/tests/ex14.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Tests all TSRK types \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown 
5*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
6c4762a1bSJed Brown   PetscInt           i, n;
7c4762a1bSJed Brown   const PetscScalar *xx;
8c4762a1bSJed Brown   /* */ PetscScalar *ff;
9c4762a1bSJed Brown 
107510d9b0SBarry Smith   PetscFunctionBeginUser;
119566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xx));
139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &ff));
14c4762a1bSJed Brown 
15*9371c9d4SSatish Balay   if (n >= 1) ff[0] = 1;
16*9371c9d4SSatish Balay   for (i = 1; i < n; i++) ff[i] = (i + 1) * (xx[i - 1] + PetscPowReal(t, i)) / 2;
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xx));
199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &ff));
20c4762a1bSJed Brown   PetscFunctionReturn(0);
21c4762a1bSJed Brown }
22c4762a1bSJed Brown 
23*9371c9d4SSatish Balay PetscErrorCode TestCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec X, PetscBool *accept) {
24c4762a1bSJed Brown   PetscInt step;
25c4762a1bSJed Brown 
267510d9b0SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
28c4762a1bSJed Brown   *accept = (step >= 2) ? PETSC_FALSE : PETSC_TRUE;
29c4762a1bSJed Brown   PetscFunctionReturn(0);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32*9371c9d4SSatish Balay static PetscErrorCode TestExplicitTS(TS ts, PetscInt order, const char subtype[]) {
33c4762a1bSJed Brown   PetscInt           i;
34c4762a1bSJed Brown   PetscReal          t;
35c4762a1bSJed Brown   Vec                U, X, Y;
36c4762a1bSJed Brown   TSType             type;
37c4762a1bSJed Brown   PetscBool          done;
38c4762a1bSJed Brown   const PetscScalar *xx;
39c4762a1bSJed Brown   const PetscScalar *yy;
40c4762a1bSJed Brown   const PetscReal    Tf  = 1;
41c4762a1bSJed Brown   const PetscReal    dt  = Tf / 8;
42c4762a1bSJed Brown   const PetscReal    eps = 100 * PETSC_MACHINE_EPSILON;
43c4762a1bSJed Brown   TSAdapt            adapt;
44c4762a1bSJed Brown   PetscInt           step;
45c4762a1bSJed Brown   TSConvergedReason  reason;
46c4762a1bSJed Brown 
477510d9b0SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(TSGetType(ts, &type));
499566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
509566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(U));
519566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
529566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0));
539566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
549566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, Tf));
559566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
569566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, NULL));
579566063dSJacob Faibussowitsch   PetscCall(TSRollBack(ts));
589566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, NULL));
599566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
629566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &X));
639566063dSJacob Faibussowitsch   PetscCall(TSEvaluateStep(ts, order, X, NULL));
649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xx));
65c4762a1bSJed Brown   for (i = 0; i < order; i++) {
66c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(xx[i]) - PetscPowReal(t, i + 1));
673c633725SBarry Smith     PetscCheck(error <= eps, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad solution, error %g, %s '%s'", (double)error, type, subtype);
68c4762a1bSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xx));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Y));
749566063dSJacob Faibussowitsch   PetscCall(TSEvaluateStep(ts, order - 1, Y, &done));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &yy));
76c4762a1bSJed Brown   for (i = 0; done && i < order - 1; i++) {
77c4762a1bSJed Brown     PetscReal error = PetscAbsReal(PetscRealPart(yy[i]) - PetscPowReal(t, i + 1));
783c633725SBarry Smith     PetscCheck(error <= eps, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad estimator, error %g, %s '%s'", (double)error, type, subtype);
79c4762a1bSJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &yy));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
849566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetCheckStage(adapt, TestCheckStage));
859566063dSJacob Faibussowitsch   PetscCall(TSSetErrorIfStepFails(ts, PETSC_FALSE));
869566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
879566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0));
889566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
899566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, NULL));
909566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetCheckStage(adapt, NULL));
919566063dSJacob Faibussowitsch   PetscCall(TSSetErrorIfStepFails(ts, PETSC_TRUE));
929566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
939566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
9463a3b9bcSJacob Faibussowitsch   PetscCheck(step == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad step number %" PetscInt_FMT ", %s '%s'", step, type, subtype);
953c633725SBarry Smith   PetscCheck(reason == TS_DIVERGED_STEP_REJECTED, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad reason %s, %s '%s'", TSConvergedReasons[reason], type, subtype);
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99*9371c9d4SSatish Balay static PetscErrorCode TestTSRK(TS ts, TSRKType type) {
100c4762a1bSJed Brown   PetscInt    order;
101c4762a1bSJed Brown   TSAdapt     adapt;
102c4762a1bSJed Brown   PetscBool   rk1, rk3, rk4;
103c4762a1bSJed Brown   TSAdaptType adapttype;
104c4762a1bSJed Brown   char        savetype[32];
105c4762a1bSJed Brown 
1067510d9b0SBarry Smith   PetscFunctionBeginUser;
1079566063dSJacob Faibussowitsch   PetscCall(TSRKSetType(ts, type));
1089566063dSJacob Faibussowitsch   PetscCall(TSRKGetType(ts, &type));
1099566063dSJacob Faibussowitsch   PetscCall(TSRKGetOrder(ts, &order));
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
1129566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetType(adapt, &adapttype));
1139566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(savetype, adapttype, sizeof(savetype)));
1149566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type, TSRK1FE, &rk1));
1159566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type, TSRK3, &rk3));
1169566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(type, TSRK4, &rk4));
1179566063dSJacob Faibussowitsch   if (rk1 || rk3 || rk4) PetscCall(TSAdaptSetType(adapt, TSADAPTNONE));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(TestExplicitTS(ts, order, type));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
1229566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt, savetype));
123c4762a1bSJed Brown   PetscFunctionReturn(0);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126*9371c9d4SSatish Balay int main(int argc, char *argv[]) {
127c4762a1bSJed Brown   TS       ts;
128c4762a1bSJed Brown   Vec      X;
129c4762a1bSJed Brown   PetscInt N = 9;
130c4762a1bSJed Brown 
131327415f7SBarry Smith   PetscFunctionBeginUser;
1329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1359566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
1369566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
1379566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &X));
1389566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
1399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1409566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK1FE));
1439566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK2A));
1449566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK3));
1459566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK3BS));
1469566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK4));
1479566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK5F));
1489566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK5DP));
1499566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK5BS));
1509566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK6VR));
1519566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK7VR));
1529566063dSJacob Faibussowitsch   PetscCall(TestTSRK(ts, TSRK8VR));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(TSRollBack(ts));
1559566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
158b122ec5aSJacob Faibussowitsch   return 0;
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /*TEST
162c4762a1bSJed Brown 
163c4762a1bSJed Brown testset:
164c4762a1bSJed Brown   output_file: output/ex14.out
165c4762a1bSJed Brown   test:
166c4762a1bSJed Brown     suffix: 0
167c4762a1bSJed Brown   test:
168c4762a1bSJed Brown     suffix: 1
169c4762a1bSJed Brown     args: -ts_adapt_type none
170c4762a1bSJed Brown   test:
171c4762a1bSJed Brown     suffix: 2
172c4762a1bSJed Brown     args: -ts_adapt_type basic
173c4762a1bSJed Brown   test:
174c4762a1bSJed Brown     suffix: 3
175c4762a1bSJed Brown     args: -ts_adapt_type dsp
176c4762a1bSJed Brown 
177c4762a1bSJed Brown TEST*/
178