1 |
heimbach |
1.1 |
C $Header: /u/gcmpack/development/heimbach/ecco_env/pkg/cal/cal_stepsforday.F,v 1.4 2001/02/02 16:57:22 heimbach Exp $ |
2 |
|
|
|
3 |
|
|
#include "CAL_CPPOPTIONS.h" |
4 |
|
|
|
5 |
|
|
subroutine cal_StepsForDay( |
6 |
|
|
I iday, |
7 |
|
|
O firststep, |
8 |
|
|
O laststep, |
9 |
|
|
O nsteps, |
10 |
|
|
I mythid |
11 |
|
|
& ) |
12 |
|
|
|
13 |
|
|
c ================================================================== |
14 |
|
|
c SUBROUTINE cal_StepsForDay |
15 |
|
|
c ================================================================== |
16 |
|
|
c |
17 |
|
|
c o Given the current day of the integration this routine returns |
18 |
|
|
c first, the last and the number of model steps the will have to |
19 |
|
|
c be performed. |
20 |
|
|
c |
21 |
|
|
c This routine also checks consistency of variables quite |
22 |
|
|
c extensively. |
23 |
|
|
c |
24 |
|
|
c started: Christian Eckert eckert@mit.edu 06-Apr-2000 |
25 |
|
|
c |
26 |
|
|
c changed: |
27 |
|
|
c |
28 |
|
|
c ================================================================== |
29 |
|
|
c SUBROUTINE cal_StepsForDay |
30 |
|
|
c ================================================================== |
31 |
|
|
|
32 |
|
|
implicit none |
33 |
|
|
|
34 |
|
|
c == global variables == |
35 |
|
|
|
36 |
|
|
#include "cal.h" |
37 |
|
|
|
38 |
|
|
c == routine arguments == |
39 |
|
|
|
40 |
|
|
integer iday |
41 |
|
|
integer firststep |
42 |
|
|
integer laststep |
43 |
|
|
integer nsteps |
44 |
|
|
integer mythid |
45 |
|
|
|
46 |
|
|
c == local variables == |
47 |
|
|
|
48 |
|
|
integer ierr |
49 |
|
|
integer mdstep |
50 |
|
|
integer numdays |
51 |
|
|
integer numsteps |
52 |
|
|
integer frac1 |
53 |
|
|
integer frac2 |
54 |
|
|
integer frac3 |
55 |
|
|
integer frac4 |
56 |
|
|
integer fullsteps |
57 |
|
|
integer firstyear |
58 |
|
|
integer firstmonth |
59 |
|
|
integer firstday |
60 |
|
|
integer lyfirst |
61 |
|
|
integer startsecs |
62 |
|
|
integer lastday |
63 |
|
|
integer endsecs |
64 |
|
|
|
65 |
|
|
c == external == |
66 |
|
|
|
67 |
|
|
external cal_IntDays |
68 |
|
|
integer cal_IntDays |
69 |
|
|
external cal_IsLeap |
70 |
|
|
integer cal_IsLeap |
71 |
|
|
|
72 |
|
|
c == end of interface == |
73 |
|
|
|
74 |
|
|
numdays = cal_IntDays( mythid ) |
75 |
|
|
lyfirst = cal_IsLeap( firstyear, mythid ) |
76 |
|
|
|
77 |
|
|
mdstep = int(modelstep) |
78 |
|
|
|
79 |
|
|
firstyear = modelstartdate(1)/10000 |
80 |
|
|
firstmonth = mod(modelstartdate(1)/100,100) |
81 |
|
|
firstday = mod(modelstartdate(1),100) |
82 |
|
|
lastday = mod(modelenddate(1),100) |
83 |
|
|
|
84 |
|
|
startsecs = (modelstartdate(2)/10000)*secondsperhour + |
85 |
|
|
& mod(modelstartdate(2)/100,100)*secondsperminute + |
86 |
|
|
& mod(modelstartdate(2),100) |
87 |
|
|
endsecs = (modelenddate(2)/10000)*secondsperhour + |
88 |
|
|
& mod(modelenddate(2)/100,100)*secondsperminute + |
89 |
|
|
& mod(modelenddate(2),100) |
90 |
|
|
|
91 |
|
|
if ( numdays .eq. 1 ) then |
92 |
|
|
if ( iday .eq. firstday ) then |
93 |
|
|
c-- Get the number of steps in the first day. |
94 |
|
|
if ( firstday .eq. lastday ) then |
95 |
|
|
firststep = 1 |
96 |
|
|
laststep = modelintsteps |
97 |
|
|
else if ( mod(firstday+1,ndaymonth(firstmonth,lyfirst)) .eq. |
98 |
|
|
& lastday ) then |
99 |
|
|
c-- This can only happen if we end at midnight of the next day. |
100 |
|
|
if ( modelenddate(2) .eq. 0 ) then |
101 |
|
|
firststep = 1 |
102 |
|
|
laststep = modelintsteps |
103 |
|
|
c-- Note: This holds only if modelenddate was determined |
104 |
|
|
c-- such that it coincides with the model's final time. |
105 |
|
|
else |
106 |
|
|
c-- We do not end at midnight of the first day of |
107 |
|
|
c-- the next month. |
108 |
|
|
ierr = 2604 |
109 |
|
|
call cal_PrintError( ierr, mythid ) |
110 |
|
|
stop ' stopped in cal_StepsForDay.' |
111 |
|
|
endif |
112 |
|
|
else |
113 |
|
|
c-- The first and the last day are inconsistent with iday. |
114 |
|
|
ierr = 2603 |
115 |
|
|
call cal_PrintError( ierr, mythid ) |
116 |
|
|
stop ' stopped in cal_StepsForDay.' |
117 |
|
|
endif |
118 |
|
|
else |
119 |
|
|
c-- The variables numdays and iday are inconsistent; |
120 |
|
|
c-- ( iday .gt. numdays ). |
121 |
|
|
ierr = 2602 |
122 |
|
|
call cal_PrintError( ierr, mythid ) |
123 |
|
|
stop ' stopped in cal_StepsForDay.' |
124 |
|
|
endif |
125 |
|
|
|
126 |
|
|
else if ( numdays .gt. 1 ) then |
127 |
|
|
c-- More than one day of integration. |
128 |
|
|
if ( iday .eq. 1 ) then |
129 |
|
|
firststep = 1 |
130 |
|
|
laststep = int((secondsperday - startsecs)/mdstep) |
131 |
|
|
else if ( ( iday .gt. 1 ) .and. |
132 |
|
|
& ( iday .lt. numdays) ) then |
133 |
|
|
c-- Somewhere between first and last month. |
134 |
|
|
c-- The first steps in iday. |
135 |
|
|
fullsteps = int((secondsperday - startsecs)/mdstep) |
136 |
|
|
numsteps = fullsteps |
137 |
|
|
c-- What's left in the first day (frac1). |
138 |
|
|
frac1 = (secondsperday - startsecs) - fullsteps*mdstep |
139 |
|
|
fullsteps = int(secondsperday/modelstep) |
140 |
|
|
c-- What's left in a complete day (frac2). |
141 |
|
|
frac2 = secondsperday - fullsteps*mdstep |
142 |
|
|
c-- What's left up to the current day (frac3). |
143 |
|
|
frac3 = frac1 + frac2*(iday - 1) |
144 |
|
|
numsteps = numsteps + (iday - 1)*fullsteps + |
145 |
|
|
& frac3/mdstep |
146 |
|
|
laststep = numsteps |
147 |
|
|
firststep = laststep - secondsperday/mdstep + 1 |
148 |
|
|
|
149 |
|
|
else if ( iday .eq. numdays ) then |
150 |
|
|
c-- The last day of integration. |
151 |
|
|
c-- The first step in iday. |
152 |
|
|
fullsteps = int((secondsperday - startsecs)/mdstep) |
153 |
|
|
numsteps = fullsteps |
154 |
|
|
c-- What's left in the first day (frac1). |
155 |
|
|
frac1 = (secondsperday - startsecs) - fullsteps*mdstep |
156 |
|
|
fullsteps = int(secondsperday/modelstep) |
157 |
|
|
c-- What's left in a complete day (frac2). |
158 |
|
|
frac2 = secondsperday - fullsteps*mdstep |
159 |
|
|
c-- What's left up to the day before the last (frac3). |
160 |
|
|
frac3 = frac1 + frac2*(iday - 2) |
161 |
|
|
numsteps = numsteps + (iday - 2)*fullsteps |
162 |
|
|
c-- The last step in iday. |
163 |
|
|
if ( modelenddate(2) .eq. 0 ) then |
164 |
|
|
c-- This can only happen if we end at midnight of the next day. |
165 |
|
|
laststep = numsteps + fullsteps |
166 |
|
|
firststep = numsteps + 1 |
167 |
|
|
c-- Note: There should be no fraction left |
168 |
|
|
c-- ( mod(frac3,mdstep) = frac3/mdstep ) if modelenddate |
169 |
|
|
c-- is based on an integral number of timesteps. |
170 |
|
|
else |
171 |
|
|
frac4 = frac3 + endsecs |
172 |
|
|
numsteps = numsteps + frac4/mdstep |
173 |
|
|
laststep = numsteps |
174 |
|
|
|
175 |
|
|
c-- Note: There should be no fraction left |
176 |
|
|
c-- ( mod(frac4,mdstep = frac4/mdstep ) if modelenddate |
177 |
|
|
c-- is based on an integral number of timesteps. |
178 |
|
|
|
179 |
|
|
firststep = laststep - endsecs/mdstep + 1 |
180 |
|
|
endif |
181 |
|
|
else |
182 |
|
|
c-- The variables iday and numdays are inconsistent. |
183 |
|
|
ierr = 2605 |
184 |
|
|
call cal_PrintError( ierr, mythid ) |
185 |
|
|
stop ' stopped in cal_DaysForMonth.' |
186 |
|
|
endif |
187 |
|
|
else |
188 |
|
|
c-- The number of days to integrate is wrong; check cal_IntDays. |
189 |
|
|
ierr = 2601 |
190 |
|
|
call cal_PrintError( ierr, mythid ) |
191 |
|
|
stop ' stopped in cal_StepsForDay.' |
192 |
|
|
endif |
193 |
|
|
|
194 |
|
|
c-- The number of days to integrate in the given month. |
195 |
|
|
nsteps = laststep - firststep + 1 |
196 |
|
|
|
197 |
|
|
return |
198 |
|
|
end |
199 |
|
|
|