1c4762a1bSJed Brown static char help[] = "Tests all TSRK types \n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown
RHSFunction(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)5*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown PetscInt i, n;
8c4762a1bSJed Brown const PetscScalar *xx;
9c4762a1bSJed Brown /* */ PetscScalar *ff;
10c4762a1bSJed Brown
117510d9b0SBarry Smith PetscFunctionBeginUser;
129566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n));
139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &xx));
149566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &ff));
15c4762a1bSJed Brown
169371c9d4SSatish Balay if (n >= 1) ff[0] = 1;
179371c9d4SSatish Balay for (i = 1; i < n; i++) ff[i] = (i + 1) * (xx[i - 1] + PetscPowReal(t, i)) / 2;
18c4762a1bSJed Brown
199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &xx));
209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &ff));
213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
22c4762a1bSJed Brown }
23c4762a1bSJed Brown
TestCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec X,PetscBool * accept)24d71ae5a4SJacob Faibussowitsch PetscErrorCode TestCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec X, PetscBool *accept)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown PetscInt step;
27c4762a1bSJed Brown
287510d9b0SBarry Smith PetscFunctionBeginUser;
299566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step));
30c4762a1bSJed Brown *accept = (step >= 2) ? PETSC_FALSE : PETSC_TRUE;
313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown
TestExplicitTS(TS ts,PetscInt order,const char subtype[])34d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestExplicitTS(TS ts, PetscInt order, const char subtype[])
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown PetscInt i;
37c4762a1bSJed Brown PetscReal t;
38c4762a1bSJed Brown Vec U, X, Y;
39c4762a1bSJed Brown TSType type;
40c4762a1bSJed Brown PetscBool done;
41c4762a1bSJed Brown const PetscScalar *xx;
42c4762a1bSJed Brown const PetscScalar *yy;
43c4762a1bSJed Brown const PetscReal Tf = 1;
44c4762a1bSJed Brown const PetscReal dt = Tf / 8;
45c4762a1bSJed Brown const PetscReal eps = 100 * PETSC_MACHINE_EPSILON;
46c4762a1bSJed Brown TSAdapt adapt;
47c4762a1bSJed Brown PetscInt step;
48c4762a1bSJed Brown TSConvergedReason reason;
49c4762a1bSJed Brown
507510d9b0SBarry Smith PetscFunctionBeginUser;
519566063dSJacob Faibussowitsch PetscCall(TSGetType(ts, &type));
529566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U));
539566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(U));
549566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0));
559566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0));
569566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
579566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, Tf));
589566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
599566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL));
609566063dSJacob Faibussowitsch PetscCall(TSRollBack(ts));
619566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL));
629566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t));
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U));
659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &X));
669566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, order, X, NULL));
679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &xx));
68c4762a1bSJed Brown for (i = 0; i < order; i++) {
69c4762a1bSJed Brown PetscReal error = PetscAbsReal(PetscRealPart(xx[i]) - PetscPowReal(t, i + 1));
703c633725SBarry Smith PetscCheck(error <= eps, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad solution, error %g, %s '%s'", (double)error, type, subtype);
71c4762a1bSJed Brown }
729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &xx));
739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
74c4762a1bSJed Brown
759566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U));
769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Y));
779566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, order - 1, Y, &done));
789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &yy));
79c4762a1bSJed Brown for (i = 0; done && i < order - 1; i++) {
80c4762a1bSJed Brown PetscReal error = PetscAbsReal(PetscRealPart(yy[i]) - PetscPowReal(t, i + 1));
813c633725SBarry Smith PetscCheck(error <= eps, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad estimator, error %g, %s '%s'", (double)error, type, subtype);
82c4762a1bSJed Brown }
839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &yy));
849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y));
85c4762a1bSJed Brown
869566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
879566063dSJacob Faibussowitsch PetscCall(TSAdaptSetCheckStage(adapt, TestCheckStage));
889566063dSJacob Faibussowitsch PetscCall(TSSetErrorIfStepFails(ts, PETSC_FALSE));
899566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0));
909566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0));
919566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
929566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, NULL));
939566063dSJacob Faibussowitsch PetscCall(TSAdaptSetCheckStage(adapt, NULL));
949566063dSJacob Faibussowitsch PetscCall(TSSetErrorIfStepFails(ts, PETSC_TRUE));
959566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step));
969566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason));
9763a3b9bcSJacob Faibussowitsch PetscCheck(step == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad step number %" PetscInt_FMT ", %s '%s'", step, type, subtype);
983c633725SBarry Smith PetscCheck(reason == TS_DIVERGED_STEP_REJECTED, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad reason %s, %s '%s'", TSConvergedReasons[reason], type, subtype);
993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown
TestTSRK(TS ts,TSRKType type)102d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestTSRK(TS ts, TSRKType type)
103d71ae5a4SJacob Faibussowitsch {
104c4762a1bSJed Brown PetscInt order;
105c4762a1bSJed Brown TSAdapt adapt;
106c4762a1bSJed Brown PetscBool rk1, rk3, rk4;
107c4762a1bSJed Brown TSAdaptType adapttype;
108c4762a1bSJed Brown char savetype[32];
109c4762a1bSJed Brown
1107510d9b0SBarry Smith PetscFunctionBeginUser;
1119566063dSJacob Faibussowitsch PetscCall(TSRKSetType(ts, type));
1129566063dSJacob Faibussowitsch PetscCall(TSRKGetType(ts, &type));
1139566063dSJacob Faibussowitsch PetscCall(TSRKGetOrder(ts, &order));
114c4762a1bSJed Brown
1159566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
1169566063dSJacob Faibussowitsch PetscCall(TSAdaptGetType(adapt, &adapttype));
1179566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(savetype, adapttype, sizeof(savetype)));
1189566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type, TSRK1FE, &rk1));
1199566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type, TSRK3, &rk3));
1209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type, TSRK4, &rk4));
1219566063dSJacob Faibussowitsch if (rk1 || rk3 || rk4) PetscCall(TSAdaptSetType(adapt, TSADAPTNONE));
122c4762a1bSJed Brown
1239566063dSJacob Faibussowitsch PetscCall(TestExplicitTS(ts, order, type));
124c4762a1bSJed Brown
1259566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt));
1269566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, savetype));
1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
128c4762a1bSJed Brown }
129c4762a1bSJed Brown
main(int argc,char * argv[])130d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
131d71ae5a4SJacob Faibussowitsch {
132c4762a1bSJed Brown TS ts;
133c4762a1bSJed Brown Vec X;
134c4762a1bSJed Brown PetscInt N = 9;
135c4762a1bSJed Brown
136327415f7SBarry Smith PetscFunctionBeginUser;
1379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
138c4762a1bSJed Brown
1399566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1409566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK));
1419566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
1429566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &X));
1439566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X));
1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1459566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
146c4762a1bSJed Brown
1479566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK1FE));
1489566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK2A));
1499566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK3));
1509566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK3BS));
1519566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK4));
1529566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK5F));
1539566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK5DP));
1549566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK5BS));
1559566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK6VR));
1569566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK7VR));
1579566063dSJacob Faibussowitsch PetscCall(TestTSRK(ts, TSRK8VR));
158c4762a1bSJed Brown
1599566063dSJacob Faibussowitsch PetscCall(TSRollBack(ts));
1609566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
161c4762a1bSJed Brown
1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
163b122ec5aSJacob Faibussowitsch return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
166c4762a1bSJed Brown /*TEST
167c4762a1bSJed Brown
168c4762a1bSJed Brown testset:
1693886731fSPierre Jolivet output_file: output/empty.out
170c4762a1bSJed Brown test:
171c4762a1bSJed Brown suffix: 0
1729d5502f9SJunchao Zhang requires: !single
173c4762a1bSJed Brown test:
174c4762a1bSJed Brown suffix: 1
175c4762a1bSJed Brown args: -ts_adapt_type none
176c4762a1bSJed Brown test:
177c4762a1bSJed Brown suffix: 2
1789d5502f9SJunchao Zhang requires: !single
179c4762a1bSJed Brown args: -ts_adapt_type basic
180c4762a1bSJed Brown test:
181c4762a1bSJed Brown suffix: 3
182c4762a1bSJed Brown args: -ts_adapt_type dsp
183c4762a1bSJed Brown
184c4762a1bSJed Brown TEST*/
185