1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2bc952696SBarry Smith
3f3334a03SStefano Zampini typedef struct {
4f3334a03SStefano Zampini PetscViewer viewer;
5f3334a03SStefano Zampini } TSTrajectory_Basic;
6533251d3SHong Zhang /*
7533251d3SHong Zhang For n-th time step, TSTrajectorySet_Basic always saves the solution X(t_n) and the current time t_n,
8533251d3SHong Zhang and optionally saves the stage values Y[] between t_{n-1} and t_n, the previous time t_{n-1}, and
9533251d3SHong Zhang forward stage sensitivities S[] = dY[]/dp.
10533251d3SHong Zhang */
TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
12d71ae5a4SJacob Faibussowitsch {
13f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic *)tj->data;
149afe7f3eSBarry Smith char filename[PETSC_MAX_PATH_LEN];
156affb6f8SHong Zhang PetscInt ns, i;
166cd0a353SHong Zhang
1722a86ba0SHong Zhang PetscFunctionBegin;
189566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof(filename), tj->dirfiletemplate, stepnum));
199566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(tjbasic->viewer, filename)); /* this triggers PetscViewer to be set up again */
209566063dSJacob Faibussowitsch PetscCall(PetscViewerSetUp(tjbasic->viewer));
219566063dSJacob Faibussowitsch PetscCall(VecView(X, tjbasic->viewer));
229566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(tjbasic->viewer, &time, 1, PETSC_REAL));
23f3334a03SStefano Zampini if (stepnum && !tj->solution_only) {
24f3334a03SStefano Zampini Vec *Y;
25f3334a03SStefano Zampini PetscReal tprev;
269566063dSJacob Faibussowitsch PetscCall(TSGetStages(ts, &ns, &Y));
27bc952696SBarry Smith for (i = 0; i < ns; i++) {
28533251d3SHong Zhang /* For stiffly accurate TS methods, the last stage Y[ns-1] is the same as the solution X, thus does not need to be saved again. */
29533251d3SHong Zhang if (ts->stifflyaccurate && i == ns - 1) continue;
309566063dSJacob Faibussowitsch PetscCall(VecView(Y[i], tjbasic->viewer));
31f3334a03SStefano Zampini }
329566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev));
339566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(tjbasic->viewer, &tprev, 1, PETSC_REAL));
34f3334a03SStefano Zampini }
356affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */
366affb6f8SHong Zhang if (ts->forward_solve) {
376affb6f8SHong Zhang Mat A, *S;
386affb6f8SHong Zhang
399566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ts, NULL, &A));
409566063dSJacob Faibussowitsch PetscCall(MatView(A, tjbasic->viewer));
416affb6f8SHong Zhang if (stepnum) {
429566063dSJacob Faibussowitsch PetscCall(TSForwardGetStages(ts, &ns, &S));
4348a46eb9SPierre Jolivet for (i = 0; i < ns; i++) PetscCall(MatView(S[i], tjbasic->viewer));
446affb6f8SHong Zhang }
456affb6f8SHong Zhang }
463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
47bc952696SBarry Smith }
48bc952696SBarry Smith
TSTrajectorySetFromOptions_Basic(TSTrajectory tj,PetscOptionItems PetscOptionsObject)49ce78bad3SBarry Smith static PetscErrorCode TSTrajectorySetFromOptions_Basic(TSTrajectory tj, PetscOptionItems PetscOptionsObject)
50d71ae5a4SJacob Faibussowitsch {
51f3334a03SStefano Zampini PetscFunctionBegin;
52d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS trajectory options for Basic type");
53d0609cedSBarry Smith PetscOptionsHeadEnd();
543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55bc952696SBarry Smith }
56bc952696SBarry Smith
TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal * t)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *t)
58d71ae5a4SJacob Faibussowitsch {
59bc952696SBarry Smith PetscViewer viewer;
609afe7f3eSBarry Smith char filename[PETSC_MAX_PATH_LEN];
61f3334a03SStefano Zampini Vec Sol;
626affb6f8SHong Zhang PetscInt ns, i;
636cd0a353SHong Zhang
6422a86ba0SHong Zhang PetscFunctionBegin;
659566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof(filename), tj->dirfiletemplate, stepnum));
669566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj), filename, FILE_MODE_READ, &viewer));
679566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
689566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &Sol));
699566063dSJacob Faibussowitsch PetscCall(VecLoad(Sol, viewer));
709566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, t, 1, NULL, PETSC_REAL));
71ac1a7491SHong Zhang if (stepnum && !tj->solution_only) {
72f3334a03SStefano Zampini Vec *Y;
73f3334a03SStefano Zampini PetscReal timepre;
749566063dSJacob Faibussowitsch PetscCall(TSGetStages(ts, &ns, &Y));
756affb6f8SHong Zhang for (i = 0; i < ns; i++) {
76533251d3SHong Zhang /* For stiffly accurate TS methods, the last stage Y[ns-1] is the same as the solution X, thus does not need to be loaded again. */
77533251d3SHong Zhang if (ts->stifflyaccurate && i == ns - 1) continue;
789566063dSJacob Faibussowitsch PetscCall(VecLoad(Y[i], viewer));
79bc952696SBarry Smith }
809566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &timepre, 1, NULL, PETSC_REAL));
811baa6e33SBarry Smith if (tj->adjoint_solve_mode) PetscCall(TSSetTimeStep(ts, -(*t) + timepre));
82f3334a03SStefano Zampini }
836affb6f8SHong Zhang /* Tangent linear sensitivities needed by second-order adjoint */
846affb6f8SHong Zhang if (ts->forward_solve) {
85cc4f23bcSHong Zhang if (!ts->stifflyaccurate) {
86cc4f23bcSHong Zhang Mat A;
879566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ts, NULL, &A));
889566063dSJacob Faibussowitsch PetscCall(MatLoad(A, viewer));
89cc4f23bcSHong Zhang }
906affb6f8SHong Zhang if (stepnum) {
91cc4f23bcSHong Zhang Mat *S;
929566063dSJacob Faibussowitsch PetscCall(TSForwardGetStages(ts, &ns, &S));
9348a46eb9SPierre Jolivet for (i = 0; i < ns; i++) PetscCall(MatLoad(S[i], viewer));
946affb6f8SHong Zhang }
956affb6f8SHong Zhang }
969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
98bc952696SBarry Smith }
99bc952696SBarry Smith
TSTrajectorySetUp_Basic(TSTrajectory tj,TS ts)100d71ae5a4SJacob Faibussowitsch PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj, TS ts)
101d71ae5a4SJacob Faibussowitsch {
102f3334a03SStefano Zampini MPI_Comm comm;
103f3334a03SStefano Zampini PetscMPIInt rank;
104f3334a03SStefano Zampini
105f3334a03SStefano Zampini PetscFunctionBegin;
1069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)tj, &comm));
1079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
108dd400576SPatrick Sanan if (rank == 0) {
1098a10d460SHong Zhang char *dir = tj->dirname;
110f3334a03SStefano Zampini PetscBool flg;
111f3334a03SStefano Zampini
1128a10d460SHong Zhang if (!dir) {
1138a10d460SHong Zhang char dtempname[16] = "TS-data-XXXXXX";
1149566063dSJacob Faibussowitsch PetscCall(PetscMkdtemp(dtempname));
1159566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(dtempname, &tj->dirname));
1168a10d460SHong Zhang } else {
1179566063dSJacob Faibussowitsch PetscCall(PetscTestDirectory(dir, 'w', &flg));
118*966bd95aSPierre Jolivet PetscCheck(!flg, comm, PETSC_ERR_SUP, "Directory %s not empty", tj->dirname);
1199566063dSJacob Faibussowitsch PetscCall(PetscTestFile(dir, 'r', &flg));
1203c633725SBarry Smith PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified path is a file - not a dir: %s", dir);
1219566063dSJacob Faibussowitsch PetscCall(PetscMkdir(dir));
122f3334a03SStefano Zampini }
1238a10d460SHong Zhang }
1249566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)tj));
1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
126f3334a03SStefano Zampini }
127f3334a03SStefano Zampini
TSTrajectoryDestroy_Basic(TSTrajectory tj)128d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj)
129d71ae5a4SJacob Faibussowitsch {
130f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic *)tj->data;
13164fc91eeSBarry Smith
13264fc91eeSBarry Smith PetscFunctionBegin;
1339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&tjbasic->viewer));
1349566063dSJacob Faibussowitsch PetscCall(PetscFree(tjbasic));
1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13664fc91eeSBarry Smith }
13764fc91eeSBarry Smith
138bc952696SBarry Smith /*MC
13978fbdcc8SBarry Smith TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
14078fbdcc8SBarry Smith
1418e5aa403SBarry Smith Saves each timestep into a separate file named TS-data-XXXXXX/TS-%06d.bin. The file name can be changed.
14278fbdcc8SBarry Smith
14378fbdcc8SBarry Smith This version saves the solutions at all the stages
14478fbdcc8SBarry Smith
14578fbdcc8SBarry Smith $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
146bc952696SBarry Smith
147bc952696SBarry Smith Level: intermediate
148bc952696SBarry Smith
1491cc06b55SBarry Smith .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectory`, `TSTrajectorySetType()`, `TSTrajectorySetDirname()`, `TSTrajectorySetFile()`,
150bcf0153eSBarry Smith `TSTrajectoryType`
151bc952696SBarry Smith M*/
TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts)152d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj, TS ts)
153d71ae5a4SJacob Faibussowitsch {
154f3334a03SStefano Zampini TSTrajectory_Basic *tjbasic;
155f3334a03SStefano Zampini
156bc952696SBarry Smith PetscFunctionBegin;
1579566063dSJacob Faibussowitsch PetscCall(PetscNew(&tjbasic));
158f3334a03SStefano Zampini
1599566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)tj), &tjbasic->viewer));
1609566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(tjbasic->viewer, PETSCVIEWERBINARY));
1619566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(tjbasic->viewer, PETSC_VIEWER_NATIVE));
1629566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(tjbasic->viewer, FILE_MODE_WRITE));
163f3334a03SStefano Zampini tj->data = tjbasic;
164f3334a03SStefano Zampini
16552fbb9a2SHong Zhang tj->ops->set = TSTrajectorySet_Basic;
16652fbb9a2SHong Zhang tj->ops->get = TSTrajectoryGet_Basic;
167f3334a03SStefano Zampini tj->ops->setup = TSTrajectorySetUp_Basic;
16864fc91eeSBarry Smith tj->ops->destroy = TSTrajectoryDestroy_Basic;
169f3334a03SStefano Zampini tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic;
1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
171bc952696SBarry Smith }
172