xref: /petsc/src/ts/tests/ex17.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
14a2752e6SStefano Zampini static char help[] = "Solves the ODE du/dt = poly(t), u(0) = 0. Tests TSResize for varying size.\n\n";
24a2752e6SStefano Zampini 
34a2752e6SStefano Zampini #include <petscts.h>
44a2752e6SStefano Zampini 
poly(PetscInt p,PetscReal t)54a2752e6SStefano Zampini PetscScalar poly(PetscInt p, PetscReal t)
64a2752e6SStefano Zampini {
74a2752e6SStefano Zampini   return p ? t * poly(p - 1, t) : 1.0;
84a2752e6SStefano Zampini }
94a2752e6SStefano Zampini 
dpoly(PetscInt p,PetscReal t)104a2752e6SStefano Zampini PetscScalar dpoly(PetscInt p, PetscReal t)
114a2752e6SStefano Zampini {
124a2752e6SStefano Zampini   return p > 0 ? (PetscReal)p * poly(p - 1, t) : 0.0;
134a2752e6SStefano Zampini }
144a2752e6SStefano Zampini 
CreateVec(PetscInt lsize,Vec * out)154a2752e6SStefano Zampini PetscErrorCode CreateVec(PetscInt lsize, Vec *out)
164a2752e6SStefano Zampini {
174a2752e6SStefano Zampini   Vec x;
184a2752e6SStefano Zampini 
194a2752e6SStefano Zampini   PetscFunctionBeginUser;
204a2752e6SStefano Zampini   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
214a2752e6SStefano Zampini   PetscCall(VecSetSizes(x, lsize, PETSC_DECIDE));
224a2752e6SStefano Zampini   PetscCall(VecSetFromOptions(x));
234a2752e6SStefano Zampini   PetscCall(VecSetUp(x));
244a2752e6SStefano Zampini   *out = x;
254a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
264a2752e6SStefano Zampini }
274a2752e6SStefano Zampini 
CreateMat(PetscInt lsize,Mat * out)284a2752e6SStefano Zampini PetscErrorCode CreateMat(PetscInt lsize, Mat *out)
294a2752e6SStefano Zampini {
304a2752e6SStefano Zampini   Mat A;
314a2752e6SStefano Zampini 
324a2752e6SStefano Zampini   PetscFunctionBeginUser;
334a2752e6SStefano Zampini   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
344a2752e6SStefano Zampini   PetscCall(MatSetSizes(A, lsize, lsize, PETSC_DECIDE, PETSC_DECIDE));
354a2752e6SStefano Zampini   PetscCall(MatSetFromOptions(A));
364a2752e6SStefano Zampini   PetscCall(MatSetUp(A));
374a2752e6SStefano Zampini   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
384a2752e6SStefano Zampini   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
394a2752e6SStefano Zampini   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
404a2752e6SStefano Zampini   PetscCall(MatShift(A, (PetscReal)1));
414a2752e6SStefano Zampini   PetscCall(MatShift(A, (PetscReal)-1));
424a2752e6SStefano Zampini   *out = A;
434a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
444a2752e6SStefano Zampini }
454a2752e6SStefano Zampini 
RHSFunction(TS ts,PetscReal t,Vec x,Vec f,PetscCtx ctx)46*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec x, Vec f, PetscCtx ctx)
474a2752e6SStefano Zampini {
484a2752e6SStefano Zampini   PetscInt *order = (PetscInt *)ctx;
494a2752e6SStefano Zampini 
504a2752e6SStefano Zampini   PetscFunctionBeginUser;
514a2752e6SStefano Zampini   PetscCall(VecSet(f, dpoly(*order, t)));
524a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
534a2752e6SStefano Zampini }
544a2752e6SStefano Zampini 
RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat B,PetscCtx ctx)55*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec x, Mat A, Mat B, PetscCtx ctx)
564a2752e6SStefano Zampini {
574a2752e6SStefano Zampini   PetscFunctionBeginUser;
584a2752e6SStefano Zampini   PetscCall(MatZeroEntries(B));
594a2752e6SStefano Zampini   if (B != A) PetscCall(MatZeroEntries(A));
604a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
614a2752e6SStefano Zampini }
624a2752e6SStefano Zampini 
Transfer(TS ts,PetscInt nv,Vec vecsin[],Vec vecsout[],PetscCtx ctx)63*2a8381b2SBarry Smith PetscErrorCode Transfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], PetscCtx ctx)
644a2752e6SStefano Zampini {
654a2752e6SStefano Zampini   PetscInt n, nnew;
664a2752e6SStefano Zampini 
674a2752e6SStefano Zampini   PetscFunctionBeginUser;
684a2752e6SStefano Zampini   PetscAssert(nv > 0, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Zero vectors");
694a2752e6SStefano Zampini   PetscCall(VecGetLocalSize(vecsin[0], &n));
704a2752e6SStefano Zampini   nnew = n == 2 ? 1 : 2;
714a2752e6SStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
724a2752e6SStefano Zampini     const PetscScalar *vals;
734a2752e6SStefano Zampini 
744a2752e6SStefano Zampini     PetscCall(CreateVec(nnew, &vecsout[i]));
754a2752e6SStefano Zampini     PetscCall(VecGetArrayRead(vecsin[i], &vals));
764a2752e6SStefano Zampini     PetscCall(VecSet(vecsout[i], vals[0]));
774a2752e6SStefano Zampini     PetscCall(VecRestoreArrayRead(vecsin[i], &vals));
784a2752e6SStefano Zampini   }
794a2752e6SStefano Zampini   Mat A;
804a2752e6SStefano Zampini   PetscCall(CreateMat(nnew, &A));
814a2752e6SStefano Zampini   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
824a2752e6SStefano Zampini   PetscCall(MatDestroy(&A));
834a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
844a2752e6SStefano Zampini }
854a2752e6SStefano Zampini 
TransferSetUp(TS ts,PetscInt step,PetscReal time,Vec sol,PetscBool * resize,PetscCtx ctx)86*2a8381b2SBarry Smith PetscErrorCode TransferSetUp(TS ts, PetscInt step, PetscReal time, Vec sol, PetscBool *resize, PetscCtx ctx)
874a2752e6SStefano Zampini {
88ecc87898SStefano Zampini   PetscBool *alreadydone = (PetscBool *)ctx;
89ecc87898SStefano Zampini 
904a2752e6SStefano Zampini   PetscFunctionBeginUser;
9157508eceSPierre Jolivet   *alreadydone = (PetscBool)!*alreadydone;
92ecc87898SStefano Zampini   *resize      = *alreadydone;
934a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
944a2752e6SStefano Zampini }
954a2752e6SStefano Zampini 
Monitor(TS ts,PetscInt n,PetscReal t,Vec x,PetscCtx ctx)96*2a8381b2SBarry Smith PetscErrorCode Monitor(TS ts, PetscInt n, PetscReal t, Vec x, PetscCtx ctx)
974a2752e6SStefano Zampini {
984a2752e6SStefano Zampini   const PetscScalar *a;
994a2752e6SStefano Zampini   PetscScalar       *store = (PetscScalar *)ctx;
1004a2752e6SStefano Zampini 
1014a2752e6SStefano Zampini   PetscFunctionBeginUser;
1024a2752e6SStefano Zampini   PetscCall(VecGetArrayRead(x, &a));
1034a2752e6SStefano Zampini   if (n < 10) store[n] = a[0];
1044a2752e6SStefano Zampini   PetscCall(VecRestoreArrayRead(x, &a));
1054a2752e6SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1064a2752e6SStefano Zampini }
1074a2752e6SStefano Zampini 
main(int argc,char ** argv)1084a2752e6SStefano Zampini int main(int argc, char **argv)
1094a2752e6SStefano Zampini {
1104a2752e6SStefano Zampini   TS          ts;
1114a2752e6SStefano Zampini   Vec         x;
1124a2752e6SStefano Zampini   Mat         A;
1134a2752e6SStefano Zampini   PetscInt    order = 2;
114ecc87898SStefano Zampini   PetscScalar results[3][10];
115ecc87898SStefano Zampini   /* I would like to use 0 here, but arkimex errors with 1.e-14 discrepancy when using TSRResize without restart on some machines (mostly arm-osx) */
116ecc87898SStefano Zampini   PetscReal tol = PETSC_SMALL;
1174a2752e6SStefano Zampini 
1184a2752e6SStefano Zampini   PetscFunctionBeginUser;
119c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1204a2752e6SStefano Zampini   PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL));
1214a2752e6SStefano Zampini   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
1224a2752e6SStefano Zampini 
123ecc87898SStefano Zampini   for (PetscInt i = 0; i < 3; i++) {
124ecc87898SStefano Zampini     PetscBool alreadydone = PETSC_TRUE;
125ecc87898SStefano Zampini 
1264a2752e6SStefano Zampini     PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1274a2752e6SStefano Zampini     PetscCall(TSSetProblemType(ts, TS_LINEAR));
1284a2752e6SStefano Zampini 
1294a2752e6SStefano Zampini     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &order));
1304a2752e6SStefano Zampini 
1314a2752e6SStefano Zampini     PetscCall(CreateMat(1, &A));
1324a2752e6SStefano Zampini     PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, NULL));
1334a2752e6SStefano Zampini     PetscCall(MatDestroy(&A));
1344a2752e6SStefano Zampini 
1354a2752e6SStefano Zampini     PetscCall(CreateVec(1, &x));
1364a2752e6SStefano Zampini     PetscCall(TSSetSolution(ts, x));
1374a2752e6SStefano Zampini     PetscCall(VecDestroy(&x));
1384a2752e6SStefano Zampini 
1394a2752e6SStefano Zampini     for (PetscInt j = 0; j < 10; j++) results[i][j] = 0;
1404a2752e6SStefano Zampini     PetscCall(TSMonitorSet(ts, Monitor, results[i], NULL));
1414a2752e6SStefano Zampini     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
142ecc87898SStefano Zampini     if (i) PetscCall(TSSetResize(ts, i == 1 ? PETSC_TRUE : PETSC_FALSE, TransferSetUp, Transfer, &alreadydone));
1434a2752e6SStefano Zampini     PetscCall(TSSetTime(ts, 0));
1444a2752e6SStefano Zampini     PetscCall(TSSetTimeStep(ts, 1. / 4.));
1454a2752e6SStefano Zampini     PetscCall(TSSetMaxSteps(ts, 10));
1464a2752e6SStefano Zampini     PetscCall(TSSetFromOptions(ts));
1474a2752e6SStefano Zampini 
1484a2752e6SStefano Zampini     PetscCall(TSSolve(ts, NULL));
1494a2752e6SStefano Zampini 
1504a2752e6SStefano Zampini     PetscCall(TSDestroy(&ts));
1514a2752e6SStefano Zampini   }
1524a2752e6SStefano Zampini 
1534a2752e6SStefano Zampini   /* Dump errors if any */
1544a2752e6SStefano Zampini   PetscBool flg = PETSC_FALSE;
1554a2752e6SStefano Zampini   for (PetscInt i = 0; i < 10; i++) {
156ecc87898SStefano Zampini     PetscReal err;
157ecc87898SStefano Zampini 
158ecc87898SStefano Zampini     err = PetscAbsScalar(results[0][i] - results[1][i]);
1594a2752e6SStefano Zampini     if (err > tol) {
160ecc87898SStefano Zampini       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error with restart for step %" PetscInt_FMT ": %g\n", i, (double)err));
161ecc87898SStefano Zampini       flg = PETSC_TRUE;
162ecc87898SStefano Zampini     }
163ecc87898SStefano Zampini     err = PetscAbsScalar(results[0][i] - results[2][i]);
164ecc87898SStefano Zampini     if (err > tol) {
165ecc87898SStefano Zampini       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error without restart for step %" PetscInt_FMT ": %g\n", i, (double)err));
1664a2752e6SStefano Zampini       flg = PETSC_TRUE;
1674a2752e6SStefano Zampini     }
1684a2752e6SStefano Zampini   }
1694a2752e6SStefano Zampini   if (flg) {
1704a2752e6SStefano Zampini     PetscCall(PetscScalarView(10, results[0], PETSC_VIEWER_STDOUT_WORLD));
1714a2752e6SStefano Zampini     PetscCall(PetscScalarView(10, results[1], PETSC_VIEWER_STDOUT_WORLD));
172ecc87898SStefano Zampini     PetscCall(PetscScalarView(10, results[2], PETSC_VIEWER_STDOUT_WORLD));
1734a2752e6SStefano Zampini   }
1744a2752e6SStefano Zampini   PetscCall(PetscFinalize());
1754a2752e6SStefano Zampini   return 0;
1764a2752e6SStefano Zampini }
1774a2752e6SStefano Zampini 
1784a2752e6SStefano Zampini /*TEST
1799d5502f9SJunchao Zhang   testset:
1809d5502f9SJunchao Zhang     # using gemv gives larger error which failed error checking
1819d5502f9SJunchao Zhang     args: -vec_mdot_use_gemv 0 -vec_maxpy_use_gemv 0
1824a2752e6SStefano Zampini 
1834a2752e6SStefano Zampini     test:
1844a2752e6SStefano Zampini       suffix: bdf
1854a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type bdf -ts_bdf_order {{2 3 4 5 6}} -order 6 -ts_adapt_type {{none basic dsp}} -ksp_type preonly -pc_type lu
1863886731fSPierre Jolivet       output_file: output/empty.out
1874a2752e6SStefano Zampini 
1884a2752e6SStefano Zampini     test:
1894a2752e6SStefano Zampini       suffix: expl
1904a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type {{euler rk ssp}} -order 6 -ts_adapt_type {{none basic dsp}}
1913886731fSPierre Jolivet       output_file: output/empty.out
1924a2752e6SStefano Zampini 
1934a2752e6SStefano Zampini     test:
1944a2752e6SStefano Zampini       suffix: impl
1954a2752e6SStefano Zampini       args: -ts_adapt_wnormtype infinity -ts_type {{rosw beuler cn alpha theta arkimex}} -order 6 -ts_adapt_type {{none basic dsp}} -ksp_type preonly -pc_type lu
1963886731fSPierre Jolivet       output_file: output/empty.out
1974a2752e6SStefano Zampini 
1984a2752e6SStefano Zampini TEST*/
199