xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1bc952696SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
3bc952696SBarry Smith 
4f3334a03SStefano Zampini typedef struct {
5f3334a03SStefano Zampini   PetscViewer viewer;
6f3334a03SStefano Zampini } TSTrajectory_Basic;
7533251d3SHong Zhang /*
8533251d3SHong Zhang   For n-th time step, TSTrajectorySet_Basic always saves the solution X(t_n) and the current time t_n,
9533251d3SHong Zhang   and optionally saves the stage values Y[] between t_{n-1} and t_n, the previous time t_{n-1}, and
10533251d3SHong Zhang   forward stage sensitivities S[] = dY[]/dp.
11533251d3SHong Zhang */
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X)
13d71ae5a4SJacob Faibussowitsch {
14f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic *)tj->data;
159afe7f3eSBarry Smith   char                filename[PETSC_MAX_PATH_LEN];
166affb6f8SHong Zhang   PetscInt            ns, i;
176cd0a353SHong Zhang 
1822a86ba0SHong Zhang   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), tj->dirfiletemplate, stepnum));
209566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(tjbasic->viewer, filename)); /* this triggers PetscViewer to be set up again */
219566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetUp(tjbasic->viewer));
229566063dSJacob Faibussowitsch   PetscCall(VecView(X, tjbasic->viewer));
239566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(tjbasic->viewer, &time, 1, PETSC_REAL));
24f3334a03SStefano Zampini   if (stepnum && !tj->solution_only) {
25f3334a03SStefano Zampini     Vec      *Y;
26f3334a03SStefano Zampini     PetscReal tprev;
279566063dSJacob Faibussowitsch     PetscCall(TSGetStages(ts, &ns, &Y));
28bc952696SBarry Smith     for (i = 0; i < ns; i++) {
29533251d3SHong 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. */
30533251d3SHong Zhang       if (ts->stifflyaccurate && i == ns - 1) continue;
319566063dSJacob Faibussowitsch       PetscCall(VecView(Y[i], tjbasic->viewer));
32f3334a03SStefano Zampini     }
339566063dSJacob Faibussowitsch     PetscCall(TSGetPrevTime(ts, &tprev));
349566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(tjbasic->viewer, &tprev, 1, PETSC_REAL));
35f3334a03SStefano Zampini   }
366affb6f8SHong Zhang   /* Tangent linear sensitivities needed by second-order adjoint */
376affb6f8SHong Zhang   if (ts->forward_solve) {
386affb6f8SHong Zhang     Mat A, *S;
396affb6f8SHong Zhang 
409566063dSJacob Faibussowitsch     PetscCall(TSForwardGetSensitivities(ts, NULL, &A));
419566063dSJacob Faibussowitsch     PetscCall(MatView(A, tjbasic->viewer));
426affb6f8SHong Zhang     if (stepnum) {
439566063dSJacob Faibussowitsch       PetscCall(TSForwardGetStages(ts, &ns, &S));
4448a46eb9SPierre Jolivet       for (i = 0; i < ns; i++) PetscCall(MatView(S[i], tjbasic->viewer));
456affb6f8SHong Zhang     }
466affb6f8SHong Zhang   }
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48bc952696SBarry Smith }
49bc952696SBarry Smith 
50d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectorySetFromOptions_Basic(TSTrajectory tj, PetscOptionItems *PetscOptionsObject)
51d71ae5a4SJacob Faibussowitsch {
52f3334a03SStefano Zampini   PetscFunctionBegin;
53d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS trajectory options for Basic type");
54d0609cedSBarry Smith   PetscOptionsHeadEnd();
553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56bc952696SBarry Smith }
57bc952696SBarry Smith 
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *t)
59d71ae5a4SJacob Faibussowitsch {
60bc952696SBarry Smith   PetscViewer viewer;
619afe7f3eSBarry Smith   char        filename[PETSC_MAX_PATH_LEN];
62f3334a03SStefano Zampini   Vec         Sol;
636affb6f8SHong Zhang   PetscInt    ns, i;
646cd0a353SHong Zhang 
6522a86ba0SHong Zhang   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), tj->dirfiletemplate, stepnum));
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj), filename, FILE_MODE_READ, &viewer));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
699566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Sol));
709566063dSJacob Faibussowitsch   PetscCall(VecLoad(Sol, viewer));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, t, 1, NULL, PETSC_REAL));
72ac1a7491SHong Zhang   if (stepnum && !tj->solution_only) {
73f3334a03SStefano Zampini     Vec      *Y;
74f3334a03SStefano Zampini     PetscReal timepre;
759566063dSJacob Faibussowitsch     PetscCall(TSGetStages(ts, &ns, &Y));
766affb6f8SHong Zhang     for (i = 0; i < ns; i++) {
77533251d3SHong 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. */
78533251d3SHong Zhang       if (ts->stifflyaccurate && i == ns - 1) continue;
799566063dSJacob Faibussowitsch       PetscCall(VecLoad(Y[i], viewer));
80bc952696SBarry Smith     }
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryRead(viewer, &timepre, 1, NULL, PETSC_REAL));
821baa6e33SBarry Smith     if (tj->adjoint_solve_mode) PetscCall(TSSetTimeStep(ts, -(*t) + timepre));
83f3334a03SStefano Zampini   }
846affb6f8SHong Zhang   /* Tangent linear sensitivities needed by second-order adjoint */
856affb6f8SHong Zhang   if (ts->forward_solve) {
86cc4f23bcSHong Zhang     if (!ts->stifflyaccurate) {
87cc4f23bcSHong Zhang       Mat A;
889566063dSJacob Faibussowitsch       PetscCall(TSForwardGetSensitivities(ts, NULL, &A));
899566063dSJacob Faibussowitsch       PetscCall(MatLoad(A, viewer));
90cc4f23bcSHong Zhang     }
916affb6f8SHong Zhang     if (stepnum) {
92cc4f23bcSHong Zhang       Mat *S;
939566063dSJacob Faibussowitsch       PetscCall(TSForwardGetStages(ts, &ns, &S));
9448a46eb9SPierre Jolivet       for (i = 0; i < ns; i++) PetscCall(MatLoad(S[i], viewer));
956affb6f8SHong Zhang     }
966affb6f8SHong Zhang   }
979566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99bc952696SBarry Smith }
100bc952696SBarry Smith 
101d71ae5a4SJacob Faibussowitsch PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj, TS ts)
102d71ae5a4SJacob Faibussowitsch {
103f3334a03SStefano Zampini   MPI_Comm    comm;
104f3334a03SStefano Zampini   PetscMPIInt rank;
105f3334a03SStefano Zampini 
106f3334a03SStefano Zampini   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tj, &comm));
1089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
109dd400576SPatrick Sanan   if (rank == 0) {
1108a10d460SHong Zhang     char     *dir = tj->dirname;
111f3334a03SStefano Zampini     PetscBool flg;
112f3334a03SStefano Zampini 
1138a10d460SHong Zhang     if (!dir) {
1148a10d460SHong Zhang       char dtempname[16] = "TS-data-XXXXXX";
1159566063dSJacob Faibussowitsch       PetscCall(PetscMkdtemp(dtempname));
1169566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(dtempname, &tj->dirname));
1178a10d460SHong Zhang     } else {
1189566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(dir, 'w', &flg));
119f3334a03SStefano Zampini       if (!flg) {
1209566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(dir, 'r', &flg));
1213c633725SBarry Smith         PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified path is a file - not a dir: %s", dir);
1229566063dSJacob Faibussowitsch         PetscCall(PetscMkdir(dir));
12398921bdaSJacob Faibussowitsch       } else SETERRQ(comm, PETSC_ERR_SUP, "Directory %s not empty", tj->dirname);
124f3334a03SStefano Zampini     }
1258a10d460SHong Zhang   }
1269566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)tj));
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128f3334a03SStefano Zampini }
129f3334a03SStefano Zampini 
130d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj)
131d71ae5a4SJacob Faibussowitsch {
132f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic *)tj->data;
13364fc91eeSBarry Smith 
13464fc91eeSBarry Smith   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&tjbasic->viewer));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree(tjbasic));
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13864fc91eeSBarry Smith }
13964fc91eeSBarry Smith 
140bc952696SBarry Smith /*MC
14178fbdcc8SBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
14278fbdcc8SBarry Smith 
1438e5aa403SBarry Smith       Saves each timestep into a separate file named TS-data-XXXXXX/TS-%06d.bin. The file name can be changed.
14478fbdcc8SBarry Smith 
14578fbdcc8SBarry Smith       This version saves the solutions at all the stages
14678fbdcc8SBarry Smith 
14778fbdcc8SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
148bc952696SBarry Smith 
149bc952696SBarry Smith   Level: intermediate
150bc952696SBarry Smith 
151*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSTrajectoryCreate()`, `TS`, `TSTrajectory`, `TSTrajectorySetType()`, `TSTrajectorySetDirname()`, `TSTrajectorySetFile()`,
152bcf0153eSBarry Smith           `TSTrajectoryType`
153bc952696SBarry Smith M*/
154d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj, TS ts)
155d71ae5a4SJacob Faibussowitsch {
156f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic;
157f3334a03SStefano Zampini 
158bc952696SBarry Smith   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(PetscNew(&tjbasic));
160f3334a03SStefano Zampini 
1619566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)tj), &tjbasic->viewer));
1629566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(tjbasic->viewer, PETSCVIEWERBINARY));
1639566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(tjbasic->viewer, PETSC_VIEWER_NATIVE));
1649566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(tjbasic->viewer, FILE_MODE_WRITE));
165f3334a03SStefano Zampini   tj->data = tjbasic;
166f3334a03SStefano Zampini 
16752fbb9a2SHong Zhang   tj->ops->set            = TSTrajectorySet_Basic;
16852fbb9a2SHong Zhang   tj->ops->get            = TSTrajectoryGet_Basic;
169f3334a03SStefano Zampini   tj->ops->setup          = TSTrajectorySetUp_Basic;
17064fc91eeSBarry Smith   tj->ops->destroy        = TSTrajectoryDestroy_Basic;
171f3334a03SStefano Zampini   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic;
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173bc952696SBarry Smith }
174