-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtime.F90
150 lines (115 loc) · 4.29 KB
/
time.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
module times
! Module for Capreole (f90)
! Author: Garrelt Mellema
! Date: 2004-05-11
! This module needs to be checked for the F compiler
! This module contains the routine which allocates the hydro
! arrays. This way they do not end up on stack
use file_admin, only: stdinput, log_unit, file_input
use precision, only: dp
use my_mpi
use scaling, only: SCTIME
use string_manipulation, only: convert_case
use astroconstants, only: YEAR
implicit none
private
real(kind=dp),public :: frametime
integer,public :: LastFrame
real(kind=dp),public :: time
real(kind=dp),public :: dt
public :: init_time
contains
subroutine init_time (restart,restartfile)
! This routine initializes all time variables
! This may be a fresh start or a restart of a saved run
logical,intent(in) :: restart
character(len=19),intent(in) :: restartfile
character(len=10) :: str_time_unit
integer :: ierror
if (restart) then
call restart_time(restartfile,time,ierror)
#ifdef MPI
! Distribute the input parameters to the other nodes
call MPI_BCAST(time,1,MPI_DOUBLE_PRECISION,0,&
MPI_COMM_NEW,ierror)
#endif
else ! Fresh start
! Initialize the time to zero
time=0.0d0
endif
! Initialize the time step to zero
dt=0.0d0
! Ask for the input if you are processor 0.
if (rank == 0) then
if (.not.file_input) then
write (*,"(//,A,/)") "----- Output -----"
write (*,"(A,$)") "1) Time between outputs (specify unit): "
endif
read (stdinput,*) frametime,str_time_unit
if (.not.file_input) write (*,"(A,$)") "2) How many output frames: "
read (stdinput,*) LastFrame
endif
! report input parameters
if (rank == 0) then
write (log_unit,"(//,A,/)") "----- Output -----"
write (log_unit,"(A,1PE10.3,A)") "1) Time between outputs: ", &
frametime,str_time_unit
write (log_unit,"(A,I4)") "2) Number of output frames: ",LastFrame
! Convert to seconds
call convert_case(str_time_unit,0) ! conversion to lower case
select case (trim(adjustl(str_time_unit)))
case ("s","sec","second","secs","seconds")
case ("years","yrs","year","yr")
frametime=frametime*YEAR
case ("myears","myrs","myear","myr")
frametime=frametime*1e6*YEAR
case ("gyears","gyrs","gyear","gyr")
frametime=frametime*1e9*YEAR
case default
write(log_unit,*) "Time unit not recognized, assuming seconds"
end select
endif
#ifdef MPI
! Distribute the input parameters to the other nodes
call MPI_BCAST(frametime,1,MPI_DOUBLE_PRECISION,0,&
MPI_COMM_NEW,ierror)
call MPI_BCAST(LastFrame,1,MPI_INTEGER,0,MPI_COMM_NEW,ierror)
#endif
! Scale time
frametime=frametime/sctime
end subroutine init_time
!========================================================================
subroutine restart_time(filename,simtime,ierror)
! This routine retrieves the time (simtime)
! from the ah3 file filename.
! Should be called from module time
character(len=19),intent(in) :: filename ! name of ah3 file
real(kind=dp),intent(out) :: simtime ! time of restart
integer,intent(out) :: ierror
! AH3D header variables
character(len=80) :: banner
integer :: nrOfDim_in ! corresponds to parameter nrOfDim (no. of dimensions)
integer :: neq_in ! corresponds to parameter neq (no. of equations)
integer :: npr_in ! corresponds to parameter npr (no. of processors)
integer :: refinementFactor ! not used
integer :: nframe ! output counter
real(kind=dp) :: gamma_in ! corresponds to parameter gamma (adiab. index)
ierror=0
! Read in header
if (rank == 0) then
open(unit=40,file=filename,form="UNFORMATTED",status="old")
read(40) banner
read(40) nrOfDim_in
read(40) neq_in
read(40) npr_in
read(40) refinementFactor
read(40) nframe
read(40) gamma_in
read(40) simtime
close(40)
endif
write(log_unit,*) "Sim time is", simtime/YEAR," years"
! Scale time
simtime=simtime/SCTIME
end subroutine restart_time
end module times