#include /*I "petscts.h" I*/ static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer) { PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscViewerCreate(comm,viewer);CHKERRQ(ierr); ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr); ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode TSTrajectorySet_Basic(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) { PetscViewer viewer; PetscInt ns,i; Vec *Y; char filename[PETSC_MAX_PATH_LEN]; PetscReal tprev; PetscErrorCode ierr; MPI_Comm comm; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); if (stepnum == 0) { PetscMPIInt rank; ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); if (!rank) { ierr = PetscRMTree("SA-data");CHKERRQ(ierr); ierr = PetscMkdir("SA-data");CHKERRQ(ierr); } ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr); ierr = VecView(X,viewer);CHKERRQ(ierr); ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-%06d.bin",stepnum);CHKERRQ(ierr); ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr); ierr = VecView(X,viewer);CHKERRQ(ierr); ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); ierr = TSGetStages(ts,&ns,&Y);CHKERRQ(ierr); for (i=0;iops->set = TSTrajectorySet_Basic; tj->ops->get = TSTrajectoryGet_Basic; PetscFunctionReturn(0); }