/[MITgcm]/MITgcm/pkg/ctrl/ctrl_summary.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_summary.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (show annotations) (download)
Thu Nov 6 22:05:08 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint57g_pre, checkpoint57s_post, checkpoint58b_post, checkpoint57b_post, checkpoint52d_pre, checkpoint57g_post, checkpoint56b_post, checkpoint57y_post, checkpoint52j_pre, checkpoint54d_post, checkpoint54e_post, checkpoint57h_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint52l_post, checkpoint52k_post, checkpoint59, checkpoint58, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint58f_post, checkpoint52f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint54f_post, checkpoint58y_post, checkpoint58t_post, checkpoint55i_post, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint57t_post, checkpoint55c_post, checkpoint52e_pre, checkpoint57v_post, checkpoint57f_post, checkpoint52e_post, checkpoint53d_post, checkpoint57a_post, checkpoint57h_pre, checkpoint52b_pre, checkpoint54b_post, checkpoint58w_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55g_post, checkpoint52b_post, checkpoint52c_post, checkpoint58o_post, checkpoint57c_post, checkpoint58p_post, checkpoint58q_post, checkpoint52f_pre, checkpoint55d_post, checkpoint58e_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint58n_post, checkpoint57e_post, checkpoint55b_post, checkpoint53a_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint55f_post, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint52d_post, checkpoint53g_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint58k_post, checkpoint52a_pre, checkpoint58v_post, checkpoint52i_post, checkpoint52h_pre, checkpoint56a_post, checkpoint58l_post, checkpoint53f_post, checkpoint57h_done, checkpoint52j_post, checkpoint57j_post, checkpoint57f_pre, checkpoint58g_post, branch-netcdf, checkpoint58x_post, checkpoint52n_post, checkpoint53b_pre, checkpoint58h_post, checkpoint56c_post, checkpoint58j_post, checkpoint57a_pre, checkpoint55a_post, checkpoint57o_post, checkpoint57k_post, checkpoint53b_post, checkpoint52a_post, checkpoint57w_post, checkpoint58i_post, ecco_c52_e35, checkpoint57x_post, checkpoint58c_post, checkpoint58u_post, checkpoint53d_pre, checkpoint58s_post, checkpoint55e_post, checkpoint54c_post, checkpoint51u_post
Branch point for: netcdf-sm0
Changes since 1.2: +4 -6 lines
o merging from ecco-branch
o cleaned some cross-dependencies and updated CPP options

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_summary.F,v 1.2 2003/06/24 16:07:07 heimbach Exp $
2
3 #include "PACKAGES_CONFIG.h"
4 #include "CTRL_CPPOPTIONS.h"
5
6 subroutine ctrl_Summary( mythid )
7
8 c ==================================================================
9 c SUBROUTINE ctrl_Summary
10 c ==================================================================
11 c
12 c o Summarize the control vector part of the ECCO release.
13 c
14 c started: Christian Eckert eckert@mit.edu 06-Mar-2000
15 c
16 c changed: Christian Eckert eckert@mit.edu
17 c
18 c ==================================================================
19 c SUBROUTINE ctrl_Summary
20 c ==================================================================
21
22 implicit none
23
24 c == global variables ==
25
26 #include "EEPARAMS.h"
27 #include "SIZE.h"
28
29 #ifdef ALLOW_CAL
30 # include "cal.h"
31 #endif
32 #include "ctrl.h"
33
34 c == routine arguments ==
35
36 integer mythid
37
38 c == local variables ==
39
40 integer bi,bj
41 integer i,k
42 integer il
43 integer timeint(4)
44 integer nwetcenter
45 integer nwetsouth
46 integer nwetwest
47
48 character*(max_len_mbuf) msgbuf
49
50 c == external ==
51
52 integer ilnblnk
53 external ilnblnk
54
55 c == end of interface ==
56
57 write(msgbuf,'(a)')
58 &' '
59 call print_message( msgbuf, standardmessageunit,
60 & SQUEEZE_RIGHT , mythid)
61 write(msgbuf,'(a)')
62 &'// ======================================================='
63 call print_message( msgbuf, standardmessageunit,
64 & SQUEEZE_RIGHT , mythid)
65 write(msgbuf,'(a)')
66 &'// ECCO control vector configuration >>> START <<<'
67 call print_message( msgbuf, standardmessageunit,
68 & SQUEEZE_RIGHT , mythid)
69 write(msgbuf,'(a)')
70 &'// ======================================================='
71 call print_message( msgbuf, standardmessageunit,
72 & SQUEEZE_RIGHT , mythid)
73 write(msgbuf,'(a)')
74 &' '
75 call print_message( msgbuf, standardmessageunit,
76 & SQUEEZE_RIGHT , mythid)
77
78 write(msgbuf,'(a)')
79 &' Total number of ocean points per tile:'
80 call print_message( msgbuf, standardmessageunit,
81 & SQUEEZE_RIGHT , mythid)
82 write(msgbuf,'(a)')
83 &' --------------------------------------'
84 call print_message( msgbuf, standardmessageunit,
85 & SQUEEZE_RIGHT , mythid)
86 write(msgbuf,'(a)')
87 &' '
88 call print_message( msgbuf, standardmessageunit,
89 & SQUEEZE_RIGHT , mythid)
90 write(msgbuf,'(a,i8)') ' snx*sny*nr = ',snx*sny*nr
91 call print_message( msgbuf, standardmessageunit,
92 & SQUEEZE_RIGHT , mythid)
93 write(msgbuf,'(a)')
94 &' '
95 call print_message( msgbuf, standardmessageunit,
96 & SQUEEZE_RIGHT , mythid)
97 write(msgbuf,'(a)')
98 &' Number of ocean points per tile:'
99 call print_message( msgbuf, standardmessageunit,
100 & SQUEEZE_RIGHT , mythid)
101 write(msgbuf,'(a)')
102 &' --------------------------------'
103 call print_message( msgbuf, standardmessageunit,
104 & SQUEEZE_RIGHT , mythid)
105 do bj = 1,nsy
106 do bi = 1,nsx
107 nwetcenter = 0
108 nwetsouth = 0
109 nwetwest = 0
110 do k = 1,nr
111 nwetcenter = nwetcenter + nwetctile(bi,bj,k)
112 nwetsouth = nwetsouth + nwetstile(bi,bj,k)
113 nwetwest = nwetwest + nwetwtile(bi,bj,k)
114 enddo
115 write(msgbuf,'(a,i5.4,i5.4,i7.6,i7.6,i7.6)')
116 & ' bi,bj,#(c/s/w):',bi,bj,nwetcenter,
117 & nwetsouth,
118 & nwetwest
119 call print_message( msgbuf, standardmessageunit,
120 & SQUEEZE_RIGHT , mythid)
121 enddo
122 enddo
123
124 #ifdef ALLOW_THETA0_CONTROL
125 write(msgbuf,'(a)')
126 &' '
127 call print_message( msgbuf, standardmessageunit,
128 & SQUEEZE_RIGHT , mythid)
129 write(msgbuf,'(a)')
130 &' Initial state temperature contribution:'
131 call print_message( msgbuf, standardmessageunit,
132 & SQUEEZE_RIGHT , mythid)
133 write(msgbuf,'(a,i5.4)')
134 &' Control variable index: ',ncvarindex(1)
135 call print_message( msgbuf, standardmessageunit,
136 & SQUEEZE_RIGHT , mythid)
137 #endif
138 #ifdef ALLOW_SALT0_CONTROL
139 write(msgbuf,'(a)')
140 &' '
141 call print_message( msgbuf, standardmessageunit,
142 & SQUEEZE_RIGHT , mythid)
143 write(msgbuf,'(a)')
144 &' Initial state salinity contribution:'
145 call print_message( msgbuf, standardmessageunit,
146 & SQUEEZE_RIGHT , mythid)
147 write(msgbuf,'(a,i5.4)')
148 &' Control variable index: ',ncvarindex(2)
149 call print_message( msgbuf, standardmessageunit,
150 & SQUEEZE_RIGHT , mythid)
151 #endif
152 #ifdef ALLOW_HFLUX_CONTROL
153 write(msgbuf,'(a)')
154 &' '
155 call print_message( msgbuf, standardmessageunit,
156 & SQUEEZE_RIGHT , mythid)
157 write(msgbuf,'(a)')
158 &' Heat flux contribution:'
159 call print_message( msgbuf, standardmessageunit,
160 & SQUEEZE_RIGHT , mythid)
161 write(msgbuf,'(a,i5.4)')
162 &' Control variable index: ',ncvarindex(3)
163 call print_message( msgbuf, standardmessageunit,
164 & SQUEEZE_RIGHT , mythid)
165
166 il = ilnblnk(xx_hflux_file)
167 call cal_TimeInterval( xx_hfluxperiod, 'secs', timeint, mythid )
168
169 write(msgbuf,'(a)')
170 &' '
171 call print_message( msgbuf, standardmessageunit,
172 & SQUEEZE_RIGHT , mythid)
173 write(msgbuf,'(a,i9.8,i7.6,1x,a,a)')
174 &' Heat flux contribution starts at: ',
175 & (xx_hfluxstartdate(i), i=1,2),
176 & dayofweek(xx_hfluxstartdate(4)),'.'
177 call print_message( msgbuf, standardmessageunit,
178 & SQUEEZE_RIGHT , mythid)
179 write(msgbuf,'(a,i9.8,i7.6)')
180 &' Heat flux contribution period is: ',
181 & (timeint(i), i=1,2)
182 call print_message( msgbuf, standardmessageunit,
183 & SQUEEZE_RIGHT , mythid)
184 write(msgbuf,'(a)')
185 &' Heat flux contribution is read from file: '
186 call print_message( msgbuf, standardmessageunit,
187 & SQUEEZE_RIGHT , mythid)
188 write(msgbuf,'(a,a,a)')
189 &' >> ',xx_hflux_file(1:il),' <<'
190 call print_message( msgbuf, standardmessageunit,
191 & SQUEEZE_RIGHT , mythid)
192 #endif
193 #ifdef ALLOW_SFLUX_CONTROL
194 write(msgbuf,'(a)')
195 &' '
196 call print_message( msgbuf, standardmessageunit,
197 & SQUEEZE_RIGHT , mythid)
198 write(msgbuf,'(a)')
199 &' Salt flux contribution:'
200 call print_message( msgbuf, standardmessageunit,
201 & SQUEEZE_RIGHT , mythid)
202 write(msgbuf,'(a,i5.4)')
203 &' Control varibale index: ',ncvarindex(4)
204 call print_message( msgbuf, standardmessageunit,
205 & SQUEEZE_RIGHT , mythid)
206
207 il = ilnblnk(xx_sflux_file)
208 call cal_TimeInterval( xx_sfluxperiod, 'secs', timeint, mythid )
209
210 write(msgbuf,'(a)')
211 &' '
212 call print_message( msgbuf, standardmessageunit,
213 & SQUEEZE_RIGHT , mythid)
214 write(msgbuf,'(a,i9.8,i7.6,1x,a,a)')
215 &' Salt flux contribution starts at: ',
216 & (xx_sfluxstartdate(i), i=1,2),
217 & dayofweek(xx_sfluxstartdate(4)),'.'
218 call print_message( msgbuf, standardmessageunit,
219 & SQUEEZE_RIGHT , mythid)
220 write(msgbuf,'(a,i9.8,i7.6)')
221 &' Salt flux contribution period is: ',
222 & (timeint(i), i=1,2)
223 call print_message( msgbuf, standardmessageunit,
224 & SQUEEZE_RIGHT , mythid)
225 write(msgbuf,'(a)')
226 &' Salt flux contribution is read from file: '
227 call print_message( msgbuf, standardmessageunit,
228 & SQUEEZE_RIGHT , mythid)
229 write(msgbuf,'(a,a,a)')
230 &' >> ',xx_sflux_file(1:il),' <<'
231 call print_message( msgbuf, standardmessageunit,
232 & SQUEEZE_RIGHT , mythid)
233 #endif
234 #ifdef ALLOW_USTRESS_CONTROL
235 write(msgbuf,'(a)')
236 &' '
237 call print_message( msgbuf, standardmessageunit,
238 & SQUEEZE_RIGHT , mythid)
239 write(msgbuf,'(a)')
240 &' Zonal wind stress contribution:'
241 call print_message( msgbuf, standardmessageunit,
242 & SQUEEZE_RIGHT , mythid)
243 write(msgbuf,'(a,i5.4)')
244 &' Control variable index: ',ncvarindex(5)
245 call print_message( msgbuf, standardmessageunit,
246 & SQUEEZE_RIGHT , mythid)
247
248 il = ilnblnk(xx_tauu_file)
249 call cal_TimeInterval( xx_tauuperiod, 'secs', timeint, mythid )
250
251 write(msgbuf,'(a)')
252 &' '
253 call print_message( msgbuf, standardmessageunit,
254 & SQUEEZE_RIGHT , mythid)
255 write(msgbuf,'(a,i9.8,i7.6,1x,a,a)')
256 &' Zonal wind stress contribution starts at: ',
257 & (xx_tauustartdate(i), i=1,2),
258 & dayofweek(xx_tauustartdate(4)),'.'
259 call print_message( msgbuf, standardmessageunit,
260 & SQUEEZE_RIGHT , mythid)
261 write(msgbuf,'(a,i9.8,i7.6)')
262 &' Zonal wind stress contribution period is: ',
263 & (timeint(i), i=1,2)
264 call print_message( msgbuf, standardmessageunit,
265 & SQUEEZE_RIGHT , mythid)
266 write(msgbuf,'(a)')
267 &' Zonal wind stress contribution is read from file: '
268 call print_message( msgbuf, standardmessageunit,
269 & SQUEEZE_RIGHT , mythid)
270 write(msgbuf,'(a,a,a)')
271 &' >> ',xx_tauu_file(1:il),' <<'
272 call print_message( msgbuf, standardmessageunit,
273 & SQUEEZE_RIGHT , mythid)
274 #endif
275 #ifdef ALLOW_VSTRESS_CONTROL
276 write(msgbuf,'(a)')
277 &' '
278 call print_message( msgbuf, standardmessageunit,
279 & SQUEEZE_RIGHT , mythid)
280 write(msgbuf,'(a)')
281 &' Meridional wind stress contribution:'
282 call print_message( msgbuf, standardmessageunit,
283 & SQUEEZE_RIGHT , mythid)
284 write(msgbuf,'(a,i5.4)')
285 &' Control variable index: ',ncvarindex(6)
286 call print_message( msgbuf, standardmessageunit,
287 & SQUEEZE_RIGHT , mythid)
288
289 il = ilnblnk(xx_tauv_file)
290 call cal_TimeInterval( xx_tauvperiod, 'secs', timeint, mythid )
291
292 write(msgbuf,'(a)')
293 &' '
294 call print_message( msgbuf, standardmessageunit,
295 & SQUEEZE_RIGHT , mythid)
296 write(msgbuf,'(a,i9.8,i7.6,1x,a,a)')
297 &' Merid. wind stress contribution starts at: ',
298 & (xx_hfluxstartdate(i), i=1,2),
299 & dayofweek(xx_hfluxstartdate(4)),'.'
300 call print_message( msgbuf, standardmessageunit,
301 & SQUEEZE_RIGHT , mythid)
302 write(msgbuf,'(a,i9.8,i7.6)')
303 &' Merid. wind stress contribution period is: ',
304 & (timeint(i), i=1,2)
305 call print_message( msgbuf, standardmessageunit,
306 & SQUEEZE_RIGHT , mythid)
307 write(msgbuf,'(a)')
308 &' Merid. wind stress contribution is read from file: '
309 call print_message( msgbuf, standardmessageunit,
310 & SQUEEZE_RIGHT , mythid)
311 write(msgbuf,'(a,a,a)')
312 &' >> ',xx_tauv_file(1:il),' <<'
313 call print_message( msgbuf, standardmessageunit,
314 & SQUEEZE_RIGHT , mythid)
315 #endif
316
317 write(msgbuf,'(a)')
318 &' '
319 call print_message( msgbuf, standardmessageunit,
320 & SQUEEZE_RIGHT , mythid)
321 write(msgbuf,'(a)')
322 &'// ======================================================='
323 call print_message( msgbuf, standardmessageunit,
324 & SQUEEZE_RIGHT , mythid)
325 write(msgbuf,'(a)')
326 &'// ECCO control vector configuration >>> END <<<'
327 call print_message( msgbuf, standardmessageunit,
328 & SQUEEZE_RIGHT , mythid)
329 write(msgbuf,'(a)')
330 &'// ======================================================='
331 call print_message( msgbuf, standardmessageunit,
332 & SQUEEZE_RIGHT , mythid)
333 write(msgbuf,'(a)')
334 &' '
335 call print_message( msgbuf, standardmessageunit,
336 & SQUEEZE_RIGHT , mythid)
337
338 return
339 end
340

  ViewVC Help
Powered by ViewVC 1.1.22