xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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 */
12560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
13bc952696SBarry Smith {
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));
446affb6f8SHong Zhang       for (i=0; i<ns; i++) {
459566063dSJacob Faibussowitsch         PetscCall(MatView(S[i],tjbasic->viewer));
466affb6f8SHong Zhang       }
476affb6f8SHong Zhang     }
486affb6f8SHong Zhang   }
49f3334a03SStefano Zampini   PetscFunctionReturn(0);
50bc952696SBarry Smith }
51bc952696SBarry Smith 
52f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
53f3334a03SStefano Zampini {
54f3334a03SStefano Zampini   PetscFunctionBegin;
55d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"TS trajectory options for Basic type");
56d0609cedSBarry Smith   PetscOptionsHeadEnd();
57bc952696SBarry Smith   PetscFunctionReturn(0);
58bc952696SBarry Smith }
59bc952696SBarry Smith 
60560360afSLisandro Dalcin static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
61bc952696SBarry Smith {
62bc952696SBarry Smith   PetscViewer    viewer;
639afe7f3eSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
64f3334a03SStefano Zampini   Vec            Sol;
656affb6f8SHong Zhang   PetscInt       ns,i;
666cd0a353SHong Zhang 
6722a86ba0SHong Zhang   PetscFunctionBegin;
689566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum));
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE));
719566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&Sol));
729566063dSJacob Faibussowitsch   PetscCall(VecLoad(Sol,viewer));
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL));
74ac1a7491SHong Zhang   if (stepnum && !tj->solution_only) {
75f3334a03SStefano Zampini     Vec       *Y;
76f3334a03SStefano Zampini     PetscReal timepre;
779566063dSJacob Faibussowitsch     PetscCall(TSGetStages(ts,&ns,&Y));
786affb6f8SHong Zhang     for (i=0; i<ns; i++) {
79533251d3SHong 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. */
80533251d3SHong Zhang       if (ts->stifflyaccurate && i == ns-1) continue;
819566063dSJacob Faibussowitsch       PetscCall(VecLoad(Y[i],viewer));
82bc952696SBarry Smith     }
839566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL));
84*1baa6e33SBarry Smith     if (tj->adjoint_solve_mode) PetscCall(TSSetTimeStep(ts,-(*t)+timepre));
85f3334a03SStefano Zampini   }
866affb6f8SHong Zhang   /* Tangent linear sensitivities needed by second-order adjoint */
876affb6f8SHong Zhang   if (ts->forward_solve) {
886affb6f8SHong Zhang 
89cc4f23bcSHong Zhang     if (!ts->stifflyaccurate) {
90cc4f23bcSHong Zhang       Mat A;
919566063dSJacob Faibussowitsch       PetscCall(TSForwardGetSensitivities(ts,NULL,&A));
929566063dSJacob Faibussowitsch       PetscCall(MatLoad(A,viewer));
93cc4f23bcSHong Zhang     }
946affb6f8SHong Zhang     if (stepnum) {
95cc4f23bcSHong Zhang       Mat *S;
969566063dSJacob Faibussowitsch       PetscCall(TSForwardGetStages(ts,&ns,&S));
976affb6f8SHong Zhang       for (i=0; i<ns; i++) {
989566063dSJacob Faibussowitsch         PetscCall(MatLoad(S[i],viewer));
996affb6f8SHong Zhang       }
1006affb6f8SHong Zhang     }
1016affb6f8SHong Zhang   }
1029566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
103bc952696SBarry Smith   PetscFunctionReturn(0);
104bc952696SBarry Smith }
105bc952696SBarry Smith 
1068a10d460SHong Zhang PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj,TS ts)
107f3334a03SStefano Zampini {
108f3334a03SStefano Zampini   MPI_Comm       comm;
109f3334a03SStefano Zampini   PetscMPIInt    rank;
110f3334a03SStefano Zampini 
111f3334a03SStefano Zampini   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tj,&comm));
1139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
114dd400576SPatrick Sanan   if (rank == 0) {
1158a10d460SHong Zhang     char      *dir = tj->dirname;
116f3334a03SStefano Zampini     PetscBool flg;
117f3334a03SStefano Zampini 
1188a10d460SHong Zhang     if (!dir) {
1198a10d460SHong Zhang       char dtempname[16] = "TS-data-XXXXXX";
1209566063dSJacob Faibussowitsch       PetscCall(PetscMkdtemp(dtempname));
1219566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(dtempname,&tj->dirname));
1228a10d460SHong Zhang     } else {
1239566063dSJacob Faibussowitsch       PetscCall(PetscTestDirectory(dir,'w',&flg));
124f3334a03SStefano Zampini       if (!flg) {
1259566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(dir,'r',&flg));
1263c633725SBarry Smith         PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_USER,"Specified path is a file - not a dir: %s",dir);
1279566063dSJacob Faibussowitsch         PetscCall(PetscMkdir(dir));
12898921bdaSJacob Faibussowitsch       } else SETERRQ(comm,PETSC_ERR_SUP,"Directory %s not empty",tj->dirname);
129f3334a03SStefano Zampini     }
1308a10d460SHong Zhang   }
1319566063dSJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)tj));
132f3334a03SStefano Zampini   PetscFunctionReturn(0);
133f3334a03SStefano Zampini }
134f3334a03SStefano Zampini 
13564fc91eeSBarry Smith static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj)
13664fc91eeSBarry Smith {
137f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data;
13864fc91eeSBarry Smith 
13964fc91eeSBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&tjbasic->viewer));
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(tjbasic));
14264fc91eeSBarry Smith   PetscFunctionReturn(0);
14364fc91eeSBarry Smith }
14464fc91eeSBarry Smith 
145bc952696SBarry Smith /*MC
14678fbdcc8SBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
14778fbdcc8SBarry Smith 
1488e5aa403SBarry Smith       Saves each timestep into a separate file named TS-data-XXXXXX/TS-%06d.bin. The file name can be changed.
14978fbdcc8SBarry Smith 
15078fbdcc8SBarry Smith       This version saves the solutions at all the stages
15178fbdcc8SBarry Smith 
15278fbdcc8SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
153bc952696SBarry Smith 
154bc952696SBarry Smith   Level: intermediate
155bc952696SBarry Smith 
156db781477SPatrick Sanan .seealso: `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`, `TSTrajectorySetDirname()`, `TSTrajectorySetFile()`
157bc952696SBarry Smith 
158bc952696SBarry Smith M*/
159c3df6c96SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts)
160bc952696SBarry Smith {
161f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic;
162f3334a03SStefano Zampini 
163bc952696SBarry Smith   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PetscNew(&tjbasic));
165f3334a03SStefano Zampini 
1669566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjbasic->viewer));
1679566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(tjbasic->viewer,PETSCVIEWERBINARY));
1689566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(tjbasic->viewer,PETSC_VIEWER_NATIVE));
1699566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(tjbasic->viewer,FILE_MODE_WRITE));
170f3334a03SStefano Zampini   tj->data = tjbasic;
171f3334a03SStefano Zampini 
17252fbb9a2SHong Zhang   tj->ops->set            = TSTrajectorySet_Basic;
17352fbb9a2SHong Zhang   tj->ops->get            = TSTrajectoryGet_Basic;
174f3334a03SStefano Zampini   tj->ops->setup          = TSTrajectorySetUp_Basic;
17564fc91eeSBarry Smith   tj->ops->destroy        = TSTrajectoryDestroy_Basic;
176f3334a03SStefano Zampini   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic;
177bc952696SBarry Smith   PetscFunctionReturn(0);
178bc952696SBarry Smith }
179