xref: /petsc/src/ts/trajectory/impls/basic/trajbasic.c (revision ac1a74915c40b0575eedca4b18743d6557ed7a4a)
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;
7bc952696SBarry Smith 
8560360afSLisandro Dalcin static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
9bc952696SBarry Smith {
10f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data;
119afe7f3eSBarry Smith   char               filename[PETSC_MAX_PATH_LEN];
12bc952696SBarry Smith   PetscErrorCode     ierr;
136cd0a353SHong Zhang 
1422a86ba0SHong Zhang   PetscFunctionBegin;
159afe7f3eSBarry Smith   ierr = PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);CHKERRQ(ierr);
16*ac1a7491SHong Zhang   ierr = PetscViewerFileSetName(tjbasic->viewer,filename);CHKERRQ(ierr); /* this triggers PetscViewer to be set up again */
17f3334a03SStefano Zampini   ierr = PetscViewerSetUp(tjbasic->viewer);CHKERRQ(ierr);
18f3334a03SStefano Zampini   ierr = VecView(X,tjbasic->viewer);CHKERRQ(ierr);
19f3334a03SStefano Zampini   ierr = PetscViewerBinaryWrite(tjbasic->viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
20f3334a03SStefano Zampini   if (stepnum && !tj->solution_only) {
21f3334a03SStefano Zampini     Vec       *Y;
22f3334a03SStefano Zampini     PetscReal tprev;
23f3334a03SStefano Zampini     PetscInt  ns,i;
24bc952696SBarry Smith 
25c679fc15SHong Zhang     ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr);
26bc952696SBarry Smith     for (i=0;i<ns;i++) {
27f3334a03SStefano Zampini       ierr = VecView(Y[i],tjbasic->viewer);CHKERRQ(ierr);
28f3334a03SStefano Zampini     }
29f3334a03SStefano Zampini     ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr);
30f3334a03SStefano Zampini     ierr = PetscViewerBinaryWrite(tjbasic->viewer,&tprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
31f3334a03SStefano Zampini   }
32f3334a03SStefano Zampini   PetscFunctionReturn(0);
33bc952696SBarry Smith }
34bc952696SBarry Smith 
35f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
36f3334a03SStefano Zampini {
37f3334a03SStefano Zampini   PetscErrorCode ierr;
38c679fc15SHong Zhang 
39f3334a03SStefano Zampini   PetscFunctionBegin;
40f3334a03SStefano Zampini   ierr = PetscOptionsHead(PetscOptionsObject,"TS trajectory options for Basic type");CHKERRQ(ierr);
41f3334a03SStefano Zampini   ierr = PetscOptionsTail();CHKERRQ(ierr);
42bc952696SBarry Smith   PetscFunctionReturn(0);
43bc952696SBarry Smith }
44bc952696SBarry Smith 
45560360afSLisandro Dalcin static PetscErrorCode TSTrajectoryGet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
46bc952696SBarry Smith {
47bc952696SBarry Smith   PetscViewer    viewer;
489afe7f3eSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
49bc952696SBarry Smith   PetscErrorCode ierr;
50f3334a03SStefano Zampini   Vec            Sol;
516cd0a353SHong Zhang 
5222a86ba0SHong Zhang   PetscFunctionBegin;
53f3334a03SStefano Zampini   ierr = PetscSNPrintf(filename,sizeof(filename),tj->dirfiletemplate,stepnum);CHKERRQ(ierr);
54bc952696SBarry Smith   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
55bc952696SBarry Smith   ierr = TSGetSolution(ts,&Sol);CHKERRQ(ierr);
56f3334a03SStefano Zampini   ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
57bc952696SBarry Smith   ierr = VecLoad(Sol,viewer);CHKERRQ(ierr);
58c679fc15SHong Zhang   ierr = PetscViewerBinaryRead(viewer,t,1,NULL,PETSC_REAL);CHKERRQ(ierr);
59*ac1a7491SHong Zhang   if (stepnum && !tj->solution_only) {
60f3334a03SStefano Zampini     Vec       *Y;
61f3334a03SStefano Zampini     PetscInt  Nr,i;
62f3334a03SStefano Zampini     PetscReal timepre;
63bc952696SBarry Smith 
64bc952696SBarry Smith     ierr = TSGetStages(ts,&Nr,&Y);CHKERRQ(ierr);
65bc952696SBarry Smith     for (i=0; i<Nr; i++) {
66bc952696SBarry Smith       ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
67bc952696SBarry Smith     }
68c679fc15SHong Zhang     ierr = PetscViewerBinaryRead(viewer,&timepre,1,NULL,PETSC_REAL);CHKERRQ(ierr);
69f3334a03SStefano Zampini     if (tj->adjoint_solve_mode) {
70acaebbb6SHong Zhang       ierr = TSSetTimeStep(ts,-(*t)+timepre);CHKERRQ(ierr);
7126656371SHong Zhang     }
72f3334a03SStefano Zampini   }
7326656371SHong Zhang   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
74bc952696SBarry Smith   PetscFunctionReturn(0);
75bc952696SBarry Smith }
76bc952696SBarry Smith 
77f3334a03SStefano Zampini static PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory tj,TS ts)
78f3334a03SStefano Zampini {
79f3334a03SStefano Zampini   MPI_Comm       comm;
80f3334a03SStefano Zampini   PetscMPIInt    rank;
81f3334a03SStefano Zampini   PetscErrorCode ierr;
82f3334a03SStefano Zampini 
83f3334a03SStefano Zampini   PetscFunctionBegin;
84f3334a03SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)tj,&comm);CHKERRQ(ierr);
85f3334a03SStefano Zampini   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
86f3334a03SStefano Zampini   if (!rank) {
87f3334a03SStefano Zampini     const char* dir = tj->dirname;
88f3334a03SStefano Zampini     PetscBool   flg;
89f3334a03SStefano Zampini 
90f3334a03SStefano Zampini     /* I don't like running PetscRMTree on a directory */
91f3334a03SStefano Zampini     ierr = PetscTestDirectory(dir,'w',&flg);CHKERRQ(ierr);
92f3334a03SStefano Zampini     if (!flg) {
93f3334a03SStefano Zampini       ierr = PetscTestFile(dir,'r',&flg);CHKERRQ(ierr);
94f3334a03SStefano Zampini       if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Specified path is a file - not a dir: %s",dir);
95f3334a03SStefano Zampini       ierr = PetscMkdir(dir);CHKERRQ(ierr);
96*ac1a7491SHong Zhang     } else {
97*ac1a7491SHong Zhang       ierr = PetscRMTree(tj->dirname);CHKERRQ(ierr); /* avoid having to delete the folder manually */
98*ac1a7491SHong Zhang     }
99f3334a03SStefano Zampini   }
100f3334a03SStefano Zampini   ierr = PetscBarrier((PetscObject)tj);CHKERRQ(ierr);
101f3334a03SStefano Zampini   PetscFunctionReturn(0);
102f3334a03SStefano Zampini }
103f3334a03SStefano Zampini 
10464fc91eeSBarry Smith static PetscErrorCode TSTrajectoryDestroy_Basic(TSTrajectory tj)
10564fc91eeSBarry Smith {
106f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic = (TSTrajectory_Basic*)tj->data;
10764fc91eeSBarry Smith   PetscErrorCode     ierr;
10864fc91eeSBarry Smith 
10964fc91eeSBarry Smith   PetscFunctionBegin;
110f3334a03SStefano Zampini   ierr = PetscViewerDestroy(&tjbasic->viewer);CHKERRQ(ierr);
111f3334a03SStefano Zampini   ierr = PetscFree(tjbasic);CHKERRQ(ierr);
11264fc91eeSBarry Smith   PetscFunctionReturn(0);
11364fc91eeSBarry Smith }
11464fc91eeSBarry Smith 
115bc952696SBarry Smith /*MC
11678fbdcc8SBarry Smith       TSTRAJECTORYBASIC - Stores each solution of the ODE/DAE in a file
11778fbdcc8SBarry Smith 
11864e38db7SHong Zhang       Saves each timestep into a seperate file named SA-data/SA-%06d.bin. The file name can be changed.
11978fbdcc8SBarry Smith 
12078fbdcc8SBarry Smith       This version saves the solutions at all the stages
12178fbdcc8SBarry Smith 
12278fbdcc8SBarry Smith       $PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m can read in files created with this format
123bc952696SBarry Smith 
124bc952696SBarry Smith   Level: intermediate
125bc952696SBarry Smith 
12664e38db7SHong Zhang .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType(), TSTrajectorySetDirname(), TSTrajectorySetFile()
127bc952696SBarry Smith 
128bc952696SBarry Smith M*/
129c3df6c96SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Basic(TSTrajectory tj,TS ts)
130bc952696SBarry Smith {
131f3334a03SStefano Zampini   TSTrajectory_Basic *tjbasic;
132f3334a03SStefano Zampini   PetscErrorCode     ierr;
133f3334a03SStefano Zampini 
134bc952696SBarry Smith   PetscFunctionBegin;
135f3334a03SStefano Zampini   ierr = PetscNew(&tjbasic);CHKERRQ(ierr);
136f3334a03SStefano Zampini 
137f3334a03SStefano Zampini   ierr = PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjbasic->viewer);CHKERRQ(ierr);
138f3334a03SStefano Zampini   ierr = PetscViewerSetType(tjbasic->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
139f3334a03SStefano Zampini   ierr = PetscViewerPushFormat(tjbasic->viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
140f3334a03SStefano Zampini   ierr = PetscViewerFileSetMode(tjbasic->viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
141f3334a03SStefano Zampini   tj->data = tjbasic;
142f3334a03SStefano Zampini 
14352fbb9a2SHong Zhang   tj->ops->set            = TSTrajectorySet_Basic;
14452fbb9a2SHong Zhang   tj->ops->get            = TSTrajectoryGet_Basic;
145f3334a03SStefano Zampini   tj->ops->setup          = TSTrajectorySetUp_Basic;
14664fc91eeSBarry Smith   tj->ops->destroy        = TSTrajectoryDestroy_Basic;
147f3334a03SStefano Zampini   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Basic;
148bc952696SBarry Smith   PetscFunctionReturn(0);
149bc952696SBarry Smith }
150