/[MITgcm]/MITgcm/pkg/grdchk/grdchk_getadxx.F
ViewVC logotype

Contents of /MITgcm/pkg/grdchk/grdchk_getadxx.F

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


Revision 1.14 - (show annotations) (download)
Wed Aug 31 00:03:44 2005 UTC (18 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57s_post, checkpoint58b_post, checkpoint57y_post, checkpoint57r_post, checkpoint58, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint57t_post, checkpoint57v_post, checkpoint57y_pre, checkpint57u_post, checkpoint57w_post, checkpoint57x_post, checkpoint58c_post
Changes since 1.13: +3 -3 lines
Adding time-dependent SST, SSS control.

1 C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_getadxx.F,v 1.13 2005/08/06 11:02:01 heimbach Exp $
2
3 #include "CTRL_CPPOPTIONS.h"
4
5
6 subroutine grdchk_getadxx(
7 I icvrec,
8 I itile,
9 I jtile,
10 I layer,
11 I itilepos,
12 I jtilepos,
13 I xx_comp,
14 I mythid
15 & )
16
17 c ==================================================================
18 c SUBROUTINE grdchk_getadxx
19 c ==================================================================
20 c
21 c o Set component a component of the control vector; xx(loc)
22 c
23 c started: Christian Eckert eckert@mit.edu 08-Mar-2000
24 c continued: heimbach@mit.edu: 13-Jun-2001
25 c
26 c ==================================================================
27 c SUBROUTINE grdchk_getadxx
28 c ==================================================================
29
30 implicit none
31
32 c == global variables ==
33
34 #include "EEPARAMS.h"
35 #include "SIZE.h"
36 #include "ctrl.h"
37 #include "optim.h"
38 #include "grdchk.h"
39
40 c == routine arguments ==
41
42 integer icvrec
43 integer jtile
44 integer itile
45 integer layer
46 integer itilepos
47 integer jtilepos
48 _RL xx_comp
49 integer mythid
50
51 #ifdef ALLOW_GRDCHK
52 c == local variables ==
53
54 integer il
55 integer dumiter
56 _RL dumtime
57 _RL dummy
58
59 logical doglobalread
60 logical ladinit
61
62 character*(80) fname
63
64 integer i,j,k
65
66 c-- == external ==
67
68 integer ilnblnk
69 external ilnblnk
70
71 c-- == end of interface ==
72
73 doglobalread = .false.
74 ladinit = .false.
75 dumiter = 0
76 dumtime = 0. _d 0
77
78 if ( grdchkvarindex .eq. 0 ) then
79 STOP 'GRDCHK INDEX 0 NOT ALLOWED'
80
81 #ifdef ALLOW_THETA0_CONTROL
82 else if ( grdchkvarindex .eq. 1 ) then
83 il=ilnblnk( xx_theta_file )
84 write(fname(1:80),'(80a)') ' '
85 write(fname(1:80),'(3a,i10.10)')
86 & yadmark, xx_theta_file(1:il),'.',optimcycle
87
88 call active_read_xyz_loc( fname, tmpfld3d, 1,
89 & doglobalread, ladinit, optimcycle,
90 & mythid, dummy)
91
92 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
93
94 #endif /* ALLOW_THETA0_CONTROL */
95
96 #ifdef ALLOW_SALT0_CONTROL
97 else if ( grdchkvarindex .eq. 2 ) then
98 il=ilnblnk( xx_salt_file )
99 write(fname(1:80),'(80a)') ' '
100 write(fname(1:80),'(3a,i10.10)')
101 & yadmark, xx_salt_file(1:il),'.',optimcycle
102
103 call active_read_xyz_loc( fname, tmpfld3d, 1,
104 & doglobalread, ladinit, optimcycle,
105 & mythid, dummy)
106
107 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
108
109 #endif /* ALLOW_SALT0_CONTROL */
110
111 #ifdef ALLOW_HFLUX_CONTROL
112 else if ( grdchkvarindex .eq. 3 ) then
113 il=ilnblnk( xx_hflux_file )
114 write(fname(1:80),'(80a)') ' '
115 write(fname(1:80),'(3a,i10.10)')
116 & yadmark, xx_hflux_file(1:il),'.',optimcycle
117
118 call active_read_xy_loc( fname, tmpfld2d, icvrec,
119 & doglobalread, ladinit, optimcycle,
120 & mythid, dummy)
121
122 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
123
124 #endif /* ALLOW_HFLUX_CONTROL */
125
126 #ifdef ALLOW_SFLUX_CONTROL
127 else if ( grdchkvarindex .eq. 4 ) then
128 il=ilnblnk( xx_sflux_file )
129 write(fname(1:80),'(80a)') ' '
130 write(fname(1:80),'(3a,i10.10)')
131 & yadmark, xx_sflux_file(1:il),'.',optimcycle
132
133 call active_read_xy_loc( fname, tmpfld2d, icvrec,
134 & doglobalread, ladinit, optimcycle,
135 & mythid, dummy)
136
137 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
138
139 #endif /* ALLOW_SFLUX_CONTROL */
140
141 #ifdef ALLOW_USTRESS_CONTROL
142 else if ( grdchkvarindex .eq. 5 ) then
143 il=ilnblnk( xx_tauu_file )
144 write(fname(1:80),'(80a)') ' '
145 write(fname(1:80),'(3a,i10.10)')
146 & yadmark, xx_tauu_file(1:il),'.',optimcycle
147
148 call active_read_xy_loc( fname, tmpfld2d, icvrec,
149 & doglobalread, ladinit, optimcycle,
150 & mythid, dummy)
151
152 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
153
154 #endif /* ALLOW_USTRESS_CONTROL */
155
156 #ifdef ALLOW_VSTRESS_CONTROL
157 else if ( grdchkvarindex .eq. 6 ) then
158 il=ilnblnk( xx_tauv_file )
159 write(fname(1:80),'(80a)') ' '
160 write(fname(1:80),'(3a,i10.10)')
161 & yadmark, xx_tauv_file(1:il),'.',optimcycle
162
163 call active_read_xy_loc( fname, tmpfld2d, icvrec,
164 & doglobalread, ladinit, optimcycle,
165 & mythid, dummy)
166
167 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
168
169 #endif /* ALLOW_VSTRESS_CONTROL */
170
171 #ifdef ALLOW_ATEMP_CONTROL
172 else if ( grdchkvarindex .eq. 7 ) then
173 il=ilnblnk( xx_atemp_file )
174 write(fname(1:80),'(80a)') ' '
175 write(fname(1:80),'(3a,i10.10)')
176 & yadmark, xx_atemp_file(1:il),'.',optimcycle
177
178 call active_read_xy_loc( fname, tmpfld2d, icvrec,
179 & doglobalread, ladinit, optimcycle,
180 & mythid, dummy)
181
182 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
183
184 #endif /* ALLOW_ATEMP_CONTROL */
185
186 #ifdef ALLOW_AQH_CONTROL
187 else if ( grdchkvarindex .eq. 8 ) then
188 il=ilnblnk( xx_aqh_file )
189 write(fname(1:80),'(80a)') ' '
190 write(fname(1:80),'(3a,i10.10)')
191 & yadmark, xx_aqh_file(1:il),'.',optimcycle
192
193 call active_read_xy_loc( fname, tmpfld2d, icvrec,
194 & doglobalread, ladinit, optimcycle,
195 & mythid, dummy)
196
197 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
198
199 #endif /* ALLOW_AQH_CONTROL */
200
201 #ifdef ALLOW_UWIND_CONTROL
202 else if ( grdchkvarindex .eq. 9 ) then
203 il=ilnblnk( xx_uwind_file )
204 write(fname(1:80),'(80a)') ' '
205 write(fname(1:80),'(3a,i10.10)')
206 & yadmark, xx_uwind_file(1:il),'.',optimcycle
207
208 call active_read_xy_loc( fname, tmpfld2d, icvrec,
209 & doglobalread, ladinit, optimcycle,
210 & mythid, dummy)
211
212 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
213
214 #endif /* ALLOW_UWIND_CONTROL */
215
216 #ifdef ALLOW_VWIND_CONTROL
217 else if ( grdchkvarindex .eq. 10 ) then
218 il=ilnblnk( xx_vwind_file )
219 write(fname(1:80),'(80a)') ' '
220 write(fname(1:80),'(3a,i10.10)')
221 & yadmark, xx_vwind_file(1:il),'.',optimcycle
222
223 call active_read_xy_loc( fname, tmpfld2d, icvrec,
224 & doglobalread, ladinit, optimcycle,
225 & mythid, dummy)
226
227 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
228
229 #endif /* ALLOW_VWIND_CONTROL */
230
231 #ifdef ALLOW_OBCSN_CONTROL
232 else if ( grdchkvarindex .eq. 11 ) then
233 il=ilnblnk( xx_obcsn_file )
234 write(fname(1:80),'(80a)') ' '
235 write(fname(1:80),'(3a,i10.10)')
236 & yadmark, xx_obcsn_file(1:il),'.',optimcycle
237
238 call active_read_xz_loc( fname, tmpfldxz, icvrec,
239 & doglobalread, ladinit, optimcycle,
240 & mythid, dummy)
241
242 xx_comp = tmpfldxz( itilepos,layer,itile,jtile )
243
244 #endif /* ALLOW_OBCSN_CONTROL */
245
246 #ifdef ALLOW_OBCSS_CONTROL
247 else if ( grdchkvarindex .eq. 12 ) then
248 il=ilnblnk( xx_obcss_file )
249 write(fname(1:80),'(80a)') ' '
250 write(fname(1:80),'(3a,i10.10)')
251 & yadmark, xx_obcss_file(1:il),'.',optimcycle
252
253 call active_read_xz_loc( fname, tmpfldxz, icvrec,
254 & doglobalread, ladinit, optimcycle,
255 & mythid, dummy)
256
257 xx_comp = tmpfldxz( itilepos,layer,itile,jtile )
258
259 #endif /* ALLOW_OBCSS_CONTROL */
260
261 #ifdef ALLOW_OBCSW_CONTROL
262 else if ( grdchkvarindex .eq. 13 ) then
263 il=ilnblnk( xx_obcsw_file )
264 write(fname(1:80),'(80a)') ' '
265 write(fname(1:80),'(3a,i10.10)')
266 & yadmark, xx_obcsw_file(1:il),'.',optimcycle
267
268 call active_read_yz_loc( fname, tmpfldyz, icvrec,
269 & doglobalread, ladinit, optimcycle,
270 & mythid, dummy)
271
272 xx_comp = tmpfldyz( jtilepos,layer,itile,jtile )
273
274 #endif /* ALLOW_OBCSW_CONTROL */
275
276 #ifdef ALLOW_OBCSE_CONTROL
277 else if ( grdchkvarindex .eq. 14 ) then
278 il=ilnblnk( xx_obcse_file )
279 write(fname(1:80),'(80a)') ' '
280 write(fname(1:80),'(3a,i10.10)')
281 & yadmark, xx_obcse_file(1:il),'.',optimcycle
282
283 call active_read_yz_loc( fname, tmpfldyz, icvrec,
284 & doglobalread, ladinit, optimcycle,
285 & mythid, dummy)
286
287 xx_comp = tmpfldyz( jtilepos,layer,itile,jtile )
288
289 #endif /* ALLOW_OBCSE_CONTROL */
290
291 #ifdef ALLOW_TR10_CONTROL
292 else if ( grdchkvarindex .eq. 17 ) then
293 il=ilnblnk( xx_tr1_file )
294 write(fname(1:80),'(80a)') ' '
295 write(fname(1:80),'(3a,i10.10)')
296 & yadmark, xx_tr1_file(1:il),'.',optimcycle
297
298 call active_read_xyz_loc( fname, tmpfld3d, 1,
299 & doglobalread, ladinit, optimcycle,
300 & mythid, dummy)
301
302 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
303
304 #endif /* ALLOW_TR10_CONTROL */
305
306 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
307 else if ( grdchkvarindex .eq. 18 ) then
308 il=ilnblnk( xx_sst_file )
309 write(fname(1:80),'(80a)') ' '
310 write(fname(1:80),'(3a,i10.10)')
311 & yadmark, xx_sst_file(1:il),'.',optimcycle
312
313 call active_read_xy_loc( fname, tmpfld2d, icvrec,
314 & doglobalread, ladinit, optimcycle,
315 & mythid, dummy)
316
317 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
318
319 #endif /* ALLOW_SST0_CONTROL */
320
321 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
322 else if ( grdchkvarindex .eq. 19 ) then
323 il=ilnblnk( xx_sss_file )
324 write(fname(1:80),'(80a)') ' '
325 write(fname(1:80),'(3a,i10.10)')
326 & yadmark, xx_sss_file(1:il),'.',optimcycle
327
328 call active_read_xy_loc( fname, tmpfld2d, icvrec,
329 & doglobalread, ladinit, optimcycle,
330 & mythid, dummy)
331
332 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
333
334 #endif /* ALLOW_SSS0_CONTROL */
335
336 #ifdef ALLOW_HFACC_CONTROL
337 else if ( grdchkvarindex .eq. 20 ) then
338 il=ilnblnk( xx_hfacc_file )
339 write(fname(1:80),'(80a)') ' '
340 write(fname(1:80),'(3a,i10.10)')
341 & yadmark, xx_hfacc_file(1:il),'.',optimcycle
342
343 #ifdef ALLOW_HFACC3D_CONTROL
344
345 call active_read_xyz_loc( fname, tmpfld3d, icvrec,
346 & doglobalread, ladinit, optimcycle,
347 & mythid, dummy)
348
349 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
350
351 #else
352
353 call active_read_xy_loc( fname, tmpfld2d, icvrec,
354 & doglobalread, ladinit, optimcycle,
355 & mythid, dummy)
356
357 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
358
359 #endif /* ALLOW_HFACC3D_CONTROL */
360 #endif /* ALLOW_HFACC_CONTROL */
361
362 #ifdef ALLOW_EFLUXY0_CONTROL
363 else if ( grdchkvarindex .eq. 21 ) then
364 il=ilnblnk( xx_efluxy_file )
365 write(fname(1:80),'(80a)') ' '
366 write(fname(1:80),'(3a,i10.10)')
367 & yadmark, xx_efluxy_file(1:il),'.',optimcycle
368
369 call active_read_xyz_loc( fname, tmpfld3d, 1,
370 & doglobalread, ladinit, optimcycle,
371 & mythid, dummy)
372
373 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
374
375 #endif /* ALLOW_EFLUXY0_CONTROL */
376
377 #ifdef ALLOW_EFLUXP0_CONTROL
378 else if ( grdchkvarindex .eq. 22 ) then
379 il=ilnblnk( xx_efluxp_file )
380 write(fname(1:80),'(80a)') ' '
381 write(fname(1:80),'(3a,i10.10)')
382 & yadmark, xx_efluxp_file(1:il),'.',optimcycle
383
384 call active_read_xyz_loc( fname, tmpfld3d, 1,
385 & doglobalread, ladinit, optimcycle,
386 & mythid, dummy)
387
388 xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile )
389
390 #endif /* ALLOW_EFLUXP0_CONTROL */
391
392 #ifdef ALLOW_PRECIP_CONTROL
393 else if ( grdchkvarindex .eq. 32 ) then
394 il=ilnblnk( xx_precip_file )
395 write(fname(1:80),'(80a)') ' '
396 write(fname(1:80),'(3a,i10.10)')
397 & yadmark, xx_precip_file(1:il),'.',optimcycle
398
399 call active_read_xy_loc( fname, tmpfld2d, icvrec,
400 & doglobalread, ladinit, optimcycle,
401 & mythid, dummy)
402
403 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
404
405 #endif /* ALLOW_PRECIP_CONTROL */
406
407 #ifdef ALLOW_SWFLUX_CONTROL
408 else if ( grdchkvarindex .eq. 33 ) then
409 il=ilnblnk( xx_swflux_file )
410 write(fname(1:80),'(80a)') ' '
411 write(fname(1:80),'(3a,i10.10)')
412 & yadmark, xx_swflux_file(1:il),'.',optimcycle
413
414 call active_read_xy_loc( fname, tmpfld2d, icvrec,
415 & doglobalread, ladinit, optimcycle,
416 & mythid, dummy)
417
418 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
419
420 #endif /* ALLOW_SWFLUX_CONTROL */
421
422 #ifdef ALLOW_SWDOWN_CONTROL
423 else if ( grdchkvarindex .eq. 34 ) then
424 il=ilnblnk( xx_swdown_file )
425 write(fname(1:80),'(80a)') ' '
426 write(fname(1:80),'(3a,i10.10)')
427 & yadmark, xx_swdown_file(1:il),'.',optimcycle
428
429 call active_read_xy_loc( fname, tmpfld2d, icvrec,
430 & doglobalread, ladinit, optimcycle,
431 & mythid, dummy)
432
433 xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile )
434
435 #endif /* ALLOW_SWDOWN_CONTROL */
436
437 else
438 ce --> this index does not exist yet.
439 endif
440
441 #endif /* ALLOW_GRDCHK */
442
443 end
444

  ViewVC Help
Powered by ViewVC 1.1.22