xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
129371c9d4SSatish Balay static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal time, Vec X) {
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));
43*48a46eb9SPierre Jolivet       for (i = 0; i < ns; i++) PetscCall(MatView(S[i], tjbasic->viewer));
446affb6f8SHong Zhang     }
456affb6f8SHong Zhang   }
46f3334a03SStefano Zampini   PetscFunctionReturn(0);
47bc952696SBarry Smith }
48bc952696SBarry Smith 
499371c9d4SSatish Balay static PetscErrorCode TSTrajectorySetFromOptions_Basic(TSTrajectory tj, PetscOptionItems *PetscOptionsObject) {
50f3334a03SStefano Zampini   PetscFunctionBegin;
51d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS trajectory options for Basic type");
52d0609cedSBarry Smith   PetscOptionsHeadEnd();
53bc952696SBarry Smith   PetscFunctionReturn(0);
54bc952696SBarry Smith }
55bc952696SBarry Smith 
569371c9d4SSatish Balay static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *t) {
57bc952696SBarry Smith   PetscViewer viewer;
589afe7f3eSBarry Smith   char        filename[PETSC_MAX_PATH_LEN];
59f3334a03SStefano Zampini   Vec         Sol;
606affb6f8SHong Zhang   PetscInt    ns, i;
616cd0a353SHong Zhang 
6222a86ba0SHong Zhang   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename, sizeof(filename), tj->dirfiletemplate, stepnum));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj), filename, FILE_MODE_READ, &viewer));
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
669566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &Sol));
679566063dSJacob Faibussowitsch   PetscCall(VecLoad(Sol, viewer));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, t, 1, NULL, PETSC_REAL));
69ac1a7491SHong Zhang   if (stepnum && !tj->solution_only) {
70f3334a03SStefano Zampini     Vec      *Y;
71f3334a03SStefano Zampini     PetscReal timepre;
729566063dSJacob Faibussowitsch     PetscCall(TSGetStages(ts, &ns, &Y));
736affb6f8SHong Zhang     for (i = 0; i < ns; i++) {
74533251d3SHong 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. */
75533251d3SHong Zhang       if (ts->stifflyaccurate && i == ns - 1) continue;
769566063dSJacob Faibussowitsch       PetscCall(VecLoad(Y[i], viewer));
77bc952696SBarry Smith     }
789566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryRead(viewer, &timepre, 1, NULL, PETSC_REAL));
791baa6e33SBarry Smith     if (tj->adjoint_solve_mode) PetscCall(TSSetTimeStep(ts, -(*t) + timepre));
80f3334a03SStefano Zampini   }
816affb6f8SHong Zhang   /* Tangent linear sensitivities needed by second-order adjoint */
826affb6f8SHong Zhang   if (ts->forward_solve) {
83cc4f23bcSHong Zhang     if (!ts->stifflyaccurate) {
84cc4f23bcSHong Zhang       Mat A;
859566063dSJacob Faibussowitsch       PetscCall(TSForwardGetSensitivities(ts, NULL, &A));
869566063dSJacob Faibussowitsch       PetscCall(MatLoad(A, viewer));
87cc4f23bcSHong Zhang     }
886affb6f8SHong Zhang     if (stepnum) {
89cc4f23bcSHong Zhang       Mat *S;
909566063dSJacob Faibussowitsch       PetscCall(TSForwardGetStages(ts, &ns, &S));
91*48a46eb9SPierre Jolivet       for (i = 0; i < ns; i++) PetscCall(MatLoad(S[i], viewer));
926affb6f8SHong Zhang     }
936affb6f8SHong Zhang   }
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
95bc952696SBarry Smith   PetscFunctionReturn(0);
96bc952696SBarry Smith }
97bc952696SBarry Smith 
989371c9d4SSatish Balay PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj, TS ts) {
99f3334a03SStefano Zampini   MPI_Comm    comm;
100f3334a03SStefano Zampini   PetscMPIInt rank;
101f3334a03SStefano Zampini 
102f3334a03SStefano Zampini   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tj, &comm));
1049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
105dd400576SPatrick Sanan   if (rank == 0) {
1068a10d460SHong Zhang     char     *dir = tj->dirname;
107f3334a03SStefano Zampini     PetscBool flg;
108f3334a03SStefano Zampini 
1098a10d460SHong Zhang     if (!dir) {
1108a10d460SHong Zhang       char dtempname[16] = "TS-data-XXXXXX";
1119566063dSJacob Faibussowitsch       PetscCall(PetscMkdtemp(dtempname));
1129566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(dtempname, &tj->dirname));
1138a10d460SHong Zhang     } else {
1149566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(dir, 'w', &flg));
115f3334a03SStefano Zampini       if (!flg) {
1169566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(dir, 'r', &flg));
1173c633725SBarry Smith         PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_USER, "Specified path is a file - not a dir: %s", dir);
1189566063dSJacob Faibussowitsch         PetscCall(PetscMkdir(dir));
11998921bdaSJacob Faibussowitsch       } else SETERRQ(comm, PETSC_ERR_SUP, "Directory %s not empty", tj->dirname);
120f3334a03SStefano Zampini     }
1218a10d460SHong Zhang   }
1229566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)tj));
123f3334a03SStefano Zampini   PetscFunctionReturn(0);
124f3334a03SStefano Zampini }
125f3334a03SStefano Zampini 
1269371c9d4SSatish Balay static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj) {
127f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic *)tj->data;
12864fc91eeSBarry Smith 
12964fc91eeSBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&tjbasic->viewer));
1319566063dSJacob Faibussowitsch   PetscCall(PetscFree(tjbasic));
13264fc91eeSBarry Smith   PetscFunctionReturn(0);
13364fc91eeSBarry Smith }
13464fc91eeSBarry Smith 
135bc952696SBarry Smith /*MC
13678fbdcc8SBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
13778fbdcc8SBarry Smith 
1388e5aa403SBarry Smith       Saves each timestep into a separate file named TS-data-XXXXXX/TS-%06d.bin. The file name can be changed.
13978fbdcc8SBarry Smith 
14078fbdcc8SBarry Smith       This version saves the solutions at all the stages
14178fbdcc8SBarry Smith 
14278fbdcc8SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
143bc952696SBarry Smith 
144bc952696SBarry Smith   Level: intermediate
145bc952696SBarry Smith 
146db781477SPatrick Sanan .seealso: `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectorySetDirname()`, `TSTrajectorySetFile()`
147bc952696SBarry Smith 
148bc952696SBarry Smith M*/
1499371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj, TS ts) {
150f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic;
151f3334a03SStefano Zampini 
152bc952696SBarry Smith   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(PetscNew(&tjbasic));
154f3334a03SStefano Zampini 
1559566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)tj), &tjbasic->viewer));
1569566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(tjbasic->viewer, PETSCVIEWERBINARY));
1579566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(tjbasic->viewer, PETSC_VIEWER_NATIVE));
1589566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(tjbasic->viewer, FILE_MODE_WRITE));
159f3334a03SStefano Zampini   tj->data = tjbasic;
160f3334a03SStefano Zampini 
16152fbb9a2SHong Zhang   tj->ops->set            = TSTrajectorySet_Basic;
16252fbb9a2SHong Zhang   tj->ops->get            = TSTrajectoryGet_Basic;
163f3334a03SStefano Zampini   tj->ops->setup          = TSTrajectorySetUp_Basic;
16464fc91eeSBarry Smith   tj->ops->destroy        = TSTrajectoryDestroy_Basic;
165f3334a03SStefano Zampini   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic;
166bc952696SBarry Smith   PetscFunctionReturn(0);
167bc952696SBarry Smith }
168