1 |
C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_getadxx.F,v 1.21 2008/01/17 20:49:27 dfer Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_CPPOPTIONS.h" |
5 |
|
6 |
|
7 |
subroutine grdchk_getadxx( |
8 |
I icvrec, |
9 |
I itile, |
10 |
I jtile, |
11 |
I layer, |
12 |
I itilepos, |
13 |
I jtilepos, |
14 |
I xx_comp, |
15 |
I mythid |
16 |
& ) |
17 |
|
18 |
c ================================================================== |
19 |
c SUBROUTINE grdchk_getadxx |
20 |
c ================================================================== |
21 |
c |
22 |
c o Set component a component of the control vector; xx(loc) |
23 |
c |
24 |
c started: Christian Eckert eckert@mit.edu 08-Mar-2000 |
25 |
c continued: heimbach@mit.edu: 13-Jun-2001 |
26 |
c |
27 |
c ================================================================== |
28 |
c SUBROUTINE grdchk_getadxx |
29 |
c ================================================================== |
30 |
|
31 |
implicit none |
32 |
|
33 |
c == global variables == |
34 |
|
35 |
#include "EEPARAMS.h" |
36 |
#include "SIZE.h" |
37 |
#include "ctrl.h" |
38 |
#include "optim.h" |
39 |
#include "grdchk.h" |
40 |
|
41 |
c == routine arguments == |
42 |
|
43 |
integer icvrec |
44 |
integer jtile |
45 |
integer itile |
46 |
integer layer |
47 |
integer itilepos |
48 |
integer jtilepos |
49 |
_RL xx_comp |
50 |
integer mythid |
51 |
|
52 |
#ifdef ALLOW_GRDCHK |
53 |
c == local variables == |
54 |
|
55 |
integer il |
56 |
integer dumiter |
57 |
_RL dumtime |
58 |
_RL dummy |
59 |
|
60 |
logical doglobalread |
61 |
logical ladinit |
62 |
|
63 |
character*(80) fname |
64 |
|
65 |
integer i,j,k |
66 |
|
67 |
c-- == external == |
68 |
|
69 |
integer ilnblnk |
70 |
external ilnblnk |
71 |
|
72 |
c-- == end of interface == |
73 |
|
74 |
doglobalread = .false. |
75 |
ladinit = .false. |
76 |
dumiter = 0 |
77 |
dumtime = 0. _d 0 |
78 |
|
79 |
if ( grdchkvarindex .eq. 0 ) then |
80 |
STOP 'GRDCHK INDEX 0 NOT ALLOWED' |
81 |
|
82 |
#ifdef ALLOW_THETA0_CONTROL |
83 |
else if ( grdchkvarindex .eq. 1 ) then |
84 |
il=ilnblnk( xx_theta_file ) |
85 |
write(fname(1:80),'(80a)') ' ' |
86 |
write(fname(1:80),'(3a,i10.10)') |
87 |
& yadmark, xx_theta_file(1:il),'.',optimcycle |
88 |
|
89 |
call active_read_xyz( fname, tmpfld3d, 1, |
90 |
& doglobalread, ladinit, optimcycle, |
91 |
& mythid, dummy) |
92 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
93 |
#endif /* ALLOW_THETA0_CONTROL */ |
94 |
|
95 |
#ifdef ALLOW_SALT0_CONTROL |
96 |
else if ( grdchkvarindex .eq. 2 ) then |
97 |
il=ilnblnk( xx_salt_file ) |
98 |
write(fname(1:80),'(80a)') ' ' |
99 |
write(fname(1:80),'(3a,i10.10)') |
100 |
& yadmark, xx_salt_file(1:il),'.',optimcycle |
101 |
call active_read_xyz( fname, tmpfld3d, 1, |
102 |
& doglobalread, ladinit, optimcycle, |
103 |
& mythid, dummy) |
104 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
105 |
|
106 |
#endif /* ALLOW_SALT0_CONTROL */ |
107 |
|
108 |
#ifdef ALLOW_HFLUX_CONTROL |
109 |
else if ( grdchkvarindex .eq. 3 ) then |
110 |
il=ilnblnk( xx_hflux_file ) |
111 |
write(fname(1:80),'(80a)') ' ' |
112 |
write(fname(1:80),'(3a,i10.10)') |
113 |
& yadmark, xx_hflux_file(1:il),'.',optimcycle |
114 |
call active_read_xy( fname, tmpfld2d, icvrec, |
115 |
& doglobalread, ladinit, optimcycle, |
116 |
& mythid, dummy) |
117 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
118 |
#endif /* ALLOW_HFLUX_CONTROL */ |
119 |
|
120 |
#ifdef ALLOW_SFLUX_CONTROL |
121 |
else if ( grdchkvarindex .eq. 4 ) then |
122 |
il=ilnblnk( xx_sflux_file ) |
123 |
write(fname(1:80),'(80a)') ' ' |
124 |
write(fname(1:80),'(3a,i10.10)') |
125 |
& yadmark, xx_sflux_file(1:il),'.',optimcycle |
126 |
call active_read_xy( fname, tmpfld2d, icvrec, |
127 |
& doglobalread, ladinit, optimcycle, |
128 |
& mythid, dummy) |
129 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
130 |
#endif /* ALLOW_SFLUX_CONTROL */ |
131 |
|
132 |
#ifdef ALLOW_USTRESS_CONTROL |
133 |
else if ( grdchkvarindex .eq. 5 ) then |
134 |
il=ilnblnk( xx_tauu_file ) |
135 |
write(fname(1:80),'(80a)') ' ' |
136 |
write(fname(1:80),'(3a,i10.10)') |
137 |
& yadmark, xx_tauu_file(1:il),'.',optimcycle |
138 |
call active_read_xy( fname, tmpfld2d, icvrec, |
139 |
& doglobalread, ladinit, optimcycle, |
140 |
& mythid, dummy) |
141 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
142 |
#endif /* ALLOW_USTRESS_CONTROL */ |
143 |
|
144 |
#ifdef ALLOW_VSTRESS_CONTROL |
145 |
else if ( grdchkvarindex .eq. 6 ) then |
146 |
il=ilnblnk( xx_tauv_file ) |
147 |
write(fname(1:80),'(80a)') ' ' |
148 |
write(fname(1:80),'(3a,i10.10)') |
149 |
& yadmark, xx_tauv_file(1:il),'.',optimcycle |
150 |
call active_read_xy( fname, tmpfld2d, icvrec, |
151 |
& doglobalread, ladinit, optimcycle, |
152 |
& mythid, dummy) |
153 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
154 |
#endif /* ALLOW_VSTRESS_CONTROL */ |
155 |
|
156 |
#ifdef ALLOW_ATEMP_CONTROL |
157 |
else if ( grdchkvarindex .eq. 7 ) then |
158 |
il=ilnblnk( xx_atemp_file ) |
159 |
write(fname(1:80),'(80a)') ' ' |
160 |
write(fname(1:80),'(3a,i10.10)') |
161 |
& yadmark, xx_atemp_file(1:il),'.',optimcycle |
162 |
call active_read_xy( fname, tmpfld2d, icvrec, |
163 |
& doglobalread, ladinit, optimcycle, |
164 |
& mythid, dummy) |
165 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
166 |
#endif /* ALLOW_ATEMP_CONTROL */ |
167 |
|
168 |
#ifdef ALLOW_AQH_CONTROL |
169 |
else if ( grdchkvarindex .eq. 8 ) then |
170 |
il=ilnblnk( xx_aqh_file ) |
171 |
write(fname(1:80),'(80a)') ' ' |
172 |
write(fname(1:80),'(3a,i10.10)') |
173 |
& yadmark, xx_aqh_file(1:il),'.',optimcycle |
174 |
call active_read_xy( fname, tmpfld2d, icvrec, |
175 |
& doglobalread, ladinit, optimcycle, |
176 |
& mythid, dummy) |
177 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
178 |
#endif /* ALLOW_AQH_CONTROL */ |
179 |
|
180 |
#ifdef ALLOW_UWIND_CONTROL |
181 |
else if ( grdchkvarindex .eq. 9 ) then |
182 |
il=ilnblnk( xx_uwind_file ) |
183 |
write(fname(1:80),'(80a)') ' ' |
184 |
write(fname(1:80),'(3a,i10.10)') |
185 |
& yadmark, xx_uwind_file(1:il),'.',optimcycle |
186 |
call active_read_xy( fname, tmpfld2d, icvrec, |
187 |
& doglobalread, ladinit, optimcycle, |
188 |
& mythid, dummy) |
189 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
190 |
#endif /* ALLOW_UWIND_CONTROL */ |
191 |
|
192 |
#ifdef ALLOW_VWIND_CONTROL |
193 |
else if ( grdchkvarindex .eq. 10 ) then |
194 |
il=ilnblnk( xx_vwind_file ) |
195 |
write(fname(1:80),'(80a)') ' ' |
196 |
write(fname(1:80),'(3a,i10.10)') |
197 |
& yadmark, xx_vwind_file(1:il),'.',optimcycle |
198 |
call active_read_xy( fname, tmpfld2d, icvrec, |
199 |
& doglobalread, ladinit, optimcycle, |
200 |
& mythid, dummy) |
201 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
202 |
#endif /* ALLOW_VWIND_CONTROL */ |
203 |
|
204 |
#ifdef ALLOW_OBCSN_CONTROL |
205 |
else if ( grdchkvarindex .eq. 11 ) then |
206 |
il=ilnblnk( xx_obcsn_file ) |
207 |
write(fname(1:80),'(80a)') ' ' |
208 |
write(fname(1:80),'(3a,i10.10)') |
209 |
& yadmark, xx_obcsn_file(1:il),'.',optimcycle |
210 |
|
211 |
call active_read_xz( fname, tmpfldxz, icvrec, |
212 |
& doglobalread, ladinit, optimcycle, |
213 |
& mythid, dummy) |
214 |
|
215 |
xx_comp = tmpfldxz( itilepos,layer,itile,jtile ) |
216 |
|
217 |
#endif /* ALLOW_OBCSN_CONTROL */ |
218 |
|
219 |
#ifdef ALLOW_OBCSS_CONTROL |
220 |
else if ( grdchkvarindex .eq. 12 ) then |
221 |
il=ilnblnk( xx_obcss_file ) |
222 |
write(fname(1:80),'(80a)') ' ' |
223 |
write(fname(1:80),'(3a,i10.10)') |
224 |
& yadmark, xx_obcss_file(1:il),'.',optimcycle |
225 |
|
226 |
call active_read_xz( fname, tmpfldxz, icvrec, |
227 |
& doglobalread, ladinit, optimcycle, |
228 |
& mythid, dummy) |
229 |
|
230 |
xx_comp = tmpfldxz( itilepos,layer,itile,jtile ) |
231 |
|
232 |
#endif /* ALLOW_OBCSS_CONTROL */ |
233 |
|
234 |
#ifdef ALLOW_OBCSW_CONTROL |
235 |
else if ( grdchkvarindex .eq. 13 ) then |
236 |
il=ilnblnk( xx_obcsw_file ) |
237 |
write(fname(1:80),'(80a)') ' ' |
238 |
write(fname(1:80),'(3a,i10.10)') |
239 |
& yadmark, xx_obcsw_file(1:il),'.',optimcycle |
240 |
|
241 |
call active_read_yz( fname, tmpfldyz, icvrec, |
242 |
& doglobalread, ladinit, optimcycle, |
243 |
& mythid, dummy) |
244 |
|
245 |
xx_comp = tmpfldyz( jtilepos,layer,itile,jtile ) |
246 |
|
247 |
#endif /* ALLOW_OBCSW_CONTROL */ |
248 |
|
249 |
#ifdef ALLOW_OBCSE_CONTROL |
250 |
else if ( grdchkvarindex .eq. 14 ) then |
251 |
il=ilnblnk( xx_obcse_file ) |
252 |
write(fname(1:80),'(80a)') ' ' |
253 |
write(fname(1:80),'(3a,i10.10)') |
254 |
& yadmark, xx_obcse_file(1:il),'.',optimcycle |
255 |
|
256 |
call active_read_yz( fname, tmpfldyz, icvrec, |
257 |
& doglobalread, ladinit, optimcycle, |
258 |
& mythid, dummy) |
259 |
|
260 |
xx_comp = tmpfldyz( jtilepos,layer,itile,jtile ) |
261 |
|
262 |
#endif /* ALLOW_OBCSE_CONTROL */ |
263 |
|
264 |
#ifdef ALLOW_DIFFKR_CONTROL |
265 |
else if ( grdchkvarindex .eq. 15 ) then |
266 |
il=ilnblnk( xx_diffkr_file ) |
267 |
write(fname(1:80),'(80a)') ' ' |
268 |
write(fname(1:80),'(3a,i10.10)') |
269 |
& yadmark, xx_diffkr_file(1:il),'.',optimcycle |
270 |
|
271 |
call active_read_xyz( fname, tmpfld3d, 1, |
272 |
& doglobalread, ladinit, optimcycle, |
273 |
& mythid, dummy) |
274 |
|
275 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
276 |
|
277 |
#endif /* ALLOW_DIFFKR_CONTROL */ |
278 |
|
279 |
#ifdef ALLOW_KAPGM_CONTROL |
280 |
else if ( grdchkvarindex .eq. 16 ) then |
281 |
il=ilnblnk( xx_kapgm_file ) |
282 |
write(fname(1:80),'(80a)') ' ' |
283 |
write(fname(1:80),'(3a,i10.10)') |
284 |
& yadmark, xx_kapgm_file(1:il),'.',optimcycle |
285 |
|
286 |
call active_read_xyz( fname, tmpfld3d, 1, |
287 |
& doglobalread, ladinit, optimcycle, |
288 |
& mythid, dummy) |
289 |
|
290 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
291 |
|
292 |
#endif /* ALLOW_KAPGM_CONTROL */ |
293 |
|
294 |
#ifdef ALLOW_KAPREDI_CONTROL |
295 |
else if ( grdchkvarindex .eq. 16 ) then |
296 |
il=ilnblnk( xx_kapredi_file ) |
297 |
write(fname(1:80),'(80a)') ' ' |
298 |
write(fname(1:80),'(3a,i10.10)') |
299 |
& yadmark, xx_kapredi_file(1:il),'.',optimcycle |
300 |
|
301 |
call active_read_xyz( fname, tmpfld3d, 1, |
302 |
& doglobalread, ladinit, optimcycle, |
303 |
& mythid, dummy) |
304 |
|
305 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
306 |
|
307 |
#endif /* ALLOW_KAPREDI_CONTROL */ |
308 |
|
309 |
#ifdef ALLOW_TR10_CONTROL |
310 |
else if ( grdchkvarindex .eq. 17 ) then |
311 |
il=ilnblnk( xx_tr1_file ) |
312 |
write(fname(1:80),'(80a)') ' ' |
313 |
write(fname(1:80),'(3a,i10.10)') |
314 |
& yadmark, xx_tr1_file(1:il),'.',optimcycle |
315 |
|
316 |
call active_read_xyz( fname, tmpfld3d, 1, |
317 |
& doglobalread, ladinit, optimcycle, |
318 |
& mythid, dummy) |
319 |
|
320 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
321 |
|
322 |
#endif /* ALLOW_TR10_CONTROL */ |
323 |
|
324 |
#if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL)) |
325 |
else if ( grdchkvarindex .eq. 18 ) then |
326 |
il=ilnblnk( xx_sst_file ) |
327 |
write(fname(1:80),'(80a)') ' ' |
328 |
write(fname(1:80),'(3a,i10.10)') |
329 |
& yadmark, xx_sst_file(1:il),'.',optimcycle |
330 |
call active_read_xy( fname, tmpfld2d, icvrec, |
331 |
& doglobalread, ladinit, optimcycle, |
332 |
& mythid, dummy) |
333 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
334 |
#endif /* ALLOW_SST0_CONTROL */ |
335 |
|
336 |
#if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL)) |
337 |
else if ( grdchkvarindex .eq. 19 ) then |
338 |
il=ilnblnk( xx_sss_file ) |
339 |
write(fname(1:80),'(80a)') ' ' |
340 |
write(fname(1:80),'(3a,i10.10)') |
341 |
& yadmark, xx_sss_file(1:il),'.',optimcycle |
342 |
call active_read_xy( fname, tmpfld2d, icvrec, |
343 |
& doglobalread, ladinit, optimcycle, |
344 |
& mythid, dummy) |
345 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
346 |
#endif /* ALLOW_SSS0_CONTROL */ |
347 |
|
348 |
#ifdef ALLOW_DEPTH_CONTROL |
349 |
else if ( grdchkvarindex .eq. 20 ) then |
350 |
il=ilnblnk( xx_depth_file ) |
351 |
write(fname(1:80),'(80a)') ' ' |
352 |
write(fname(1:80),'(3a,i10.10)') |
353 |
& yadmark, xx_depth_file(1:il),'.',optimcycle |
354 |
call active_read_xy( fname, tmpfld2d, icvrec, |
355 |
& doglobalread, ladinit, optimcycle, |
356 |
& mythid, dummy) |
357 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
358 |
#endif /* ALLOW_DEPTH_CONTROL */ |
359 |
|
360 |
#ifdef ALLOW_EFLUXY0_CONTROL |
361 |
else if ( grdchkvarindex .eq. 21 ) then |
362 |
il=ilnblnk( xx_efluxy_file ) |
363 |
write(fname(1:80),'(80a)') ' ' |
364 |
write(fname(1:80),'(3a,i10.10)') |
365 |
& yadmark, xx_efluxy_file(1:il),'.',optimcycle |
366 |
|
367 |
call active_read_xyz( fname, tmpfld3d, 1, |
368 |
& doglobalread, ladinit, optimcycle, |
369 |
& mythid, dummy) |
370 |
|
371 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
372 |
|
373 |
#endif /* ALLOW_EFLUXY0_CONTROL */ |
374 |
|
375 |
#ifdef ALLOW_EFLUXP0_CONTROL |
376 |
else if ( grdchkvarindex .eq. 22 ) then |
377 |
il=ilnblnk( xx_efluxp_file ) |
378 |
write(fname(1:80),'(80a)') ' ' |
379 |
write(fname(1:80),'(3a,i10.10)') |
380 |
& yadmark, xx_efluxp_file(1:il),'.',optimcycle |
381 |
|
382 |
call active_read_xyz( fname, tmpfld3d, 1, |
383 |
& doglobalread, ladinit, optimcycle, |
384 |
& mythid, dummy) |
385 |
|
386 |
xx_comp = tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) |
387 |
|
388 |
#endif /* ALLOW_EFLUXP0_CONTROL */ |
389 |
|
390 |
#ifdef ALLOW_HFLUXM_CONTROL |
391 |
else if ( grdchkvarindex .eq. 24 ) then |
392 |
il=ilnblnk( xx_hfluxm_file ) |
393 |
write(fname(1:80),'(80a)') ' ' |
394 |
write(fname(1:80),'(3a,i10.10)') |
395 |
& yadmark, xx_hfluxm_file(1:il),'.',optimcycle |
396 |
call active_read_xy( fname, tmpfld2d, icvrec, |
397 |
& doglobalread, ladinit, optimcycle, |
398 |
& mythid, dummy) |
399 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
400 |
#endif /* ALLOW_HFLUXM_CONTROL */ |
401 |
|
402 |
#ifdef ALLOW_PRECIP_CONTROL |
403 |
else if ( grdchkvarindex .eq. 32 ) then |
404 |
il=ilnblnk( xx_precip_file ) |
405 |
write(fname(1:80),'(80a)') ' ' |
406 |
write(fname(1:80),'(3a,i10.10)') |
407 |
& yadmark, xx_precip_file(1:il),'.',optimcycle |
408 |
call active_read_xy( fname, tmpfld2d, icvrec, |
409 |
& doglobalread, ladinit, optimcycle, |
410 |
& mythid, dummy) |
411 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
412 |
#endif /* ALLOW_PRECIP_CONTROL */ |
413 |
|
414 |
#ifdef ALLOW_SWFLUX_CONTROL |
415 |
else if ( grdchkvarindex .eq. 33 ) then |
416 |
il=ilnblnk( xx_swflux_file ) |
417 |
write(fname(1:80),'(80a)') ' ' |
418 |
write(fname(1:80),'(3a,i10.10)') |
419 |
& yadmark, xx_swflux_file(1:il),'.',optimcycle |
420 |
call active_read_xy( fname, tmpfld2d, icvrec, |
421 |
& doglobalread, ladinit, optimcycle, |
422 |
& mythid, dummy) |
423 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
424 |
#endif /* ALLOW_SWFLUX_CONTROL */ |
425 |
|
426 |
#ifdef ALLOW_SWDOWN_CONTROL |
427 |
else if ( grdchkvarindex .eq. 34 ) then |
428 |
il=ilnblnk( xx_swdown_file ) |
429 |
write(fname(1:80),'(80a)') ' ' |
430 |
write(fname(1:80),'(3a,i10.10)') |
431 |
& yadmark, xx_swdown_file(1:il),'.',optimcycle |
432 |
call active_read_xy( fname, tmpfld2d, icvrec, |
433 |
& doglobalread, ladinit, optimcycle, |
434 |
& mythid, dummy) |
435 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
436 |
#endif /* ALLOW_SWDOWN_CONTROL */ |
437 |
|
438 |
#ifdef ALLOW_LWFLUX_CONTROL |
439 |
else if ( grdchkvarindex .eq. 35 ) then |
440 |
il=ilnblnk( xx_lwflux_file ) |
441 |
write(fname(1:80),'(80a)') ' ' |
442 |
write(fname(1:80),'(3a,i10.10)') |
443 |
& yadmark, xx_lwflux_file(1:il),'.',optimcycle |
444 |
call active_read_xy( fname, tmpfld2d, icvrec, |
445 |
& doglobalread, ladinit, optimcycle, |
446 |
& mythid, dummy) |
447 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
448 |
#endif /* ALLOW_LWFLUX_CONTROL */ |
449 |
|
450 |
#ifdef ALLOW_LWDOWN_CONTROL |
451 |
else if ( grdchkvarindex .eq. 36 ) then |
452 |
il=ilnblnk( xx_lwdown_file ) |
453 |
write(fname(1:80),'(80a)') ' ' |
454 |
write(fname(1:80),'(3a,i10.10)') |
455 |
& yadmark, xx_lwdown_file(1:il),'.',optimcycle |
456 |
call active_read_xy( fname, tmpfld2d, icvrec, |
457 |
& doglobalread, ladinit, optimcycle, |
458 |
& mythid, dummy) |
459 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
460 |
#endif /* ALLOW_LWDOWN_CONTROL */ |
461 |
|
462 |
#ifdef ALLOW_EVAP_CONTROL |
463 |
else if ( grdchkvarindex .eq. 37 ) then |
464 |
il=ilnblnk( xx_evap_file ) |
465 |
write(fname(1:80),'(80a)') ' ' |
466 |
write(fname(1:80),'(3a,i10.10)') |
467 |
& yadmark, xx_evap_file(1:il),'.',optimcycle |
468 |
call active_read_xy( fname, tmpfld2d, icvrec, |
469 |
& doglobalread, ladinit, optimcycle, |
470 |
& mythid, dummy) |
471 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
472 |
#endif /* ALLOW_EVAP_CONTROL */ |
473 |
|
474 |
#ifdef ALLOW_SNOWPRECIP_CONTROL |
475 |
else if ( grdchkvarindex .eq. 38 ) then |
476 |
il=ilnblnk( xx_snowprecip_file ) |
477 |
write(fname(1:80),'(80a)') ' ' |
478 |
write(fname(1:80),'(3a,i10.10)') |
479 |
& yadmark, xx_snowprecip_file(1:il),'.',optimcycle |
480 |
call active_read_xy( fname, tmpfld2d, icvrec, |
481 |
& doglobalread, ladinit, optimcycle, |
482 |
& mythid, dummy) |
483 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
484 |
#endif /* ALLOW_SNOWPRECIP_CONTROL */ |
485 |
|
486 |
#ifdef ALLOW_APRESSURE_CONTROL |
487 |
else if ( grdchkvarindex .eq. 39 ) then |
488 |
il=ilnblnk( xx_apressure_file ) |
489 |
write(fname(1:80),'(80a)') ' ' |
490 |
write(fname(1:80),'(3a,i10.10)') |
491 |
& yadmark, xx_apressure_file(1:il),'.',optimcycle |
492 |
|
493 |
call active_read_xy( fname, tmpfld2d, icvrec, |
494 |
& doglobalread, ladinit, optimcycle, |
495 |
& mythid, dummy) |
496 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
497 |
#endif /* ALLOW_APRESSURE_CONTROL */ |
498 |
|
499 |
#ifdef ALLOW_RUNOFF_CONTROL |
500 |
else if ( grdchkvarindex .eq. 40 ) then |
501 |
il=ilnblnk( xx_runoff_file ) |
502 |
write(fname(1:80),'(80a)') ' ' |
503 |
write(fname(1:80),'(3a,i10.10)') |
504 |
& yadmark, xx_runoff_file(1:il),'.',optimcycle |
505 |
call active_read_xy( fname, tmpfld2d, icvrec, |
506 |
& doglobalread, ladinit, optimcycle, |
507 |
& mythid, dummy) |
508 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
509 |
#endif /* ALLOW_RUNOFF_CONTROL */ |
510 |
|
511 |
#ifdef ALLOW_SIAREA_CONTROL |
512 |
else if ( grdchkvarindex .eq. 41 ) then |
513 |
il=ilnblnk( xx_siarea_file ) |
514 |
write(fname(1:80),'(80a)') ' ' |
515 |
write(fname(1:80),'(3a,i10.10)') |
516 |
& yadmark, xx_siarea_file(1:il),'.',optimcycle |
517 |
call active_read_xy( fname, tmpfld2d, icvrec, |
518 |
& doglobalread, ladinit, optimcycle, |
519 |
& mythid, dummy) |
520 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
521 |
#endif /* ALLOW_SIAREA_CONTROL */ |
522 |
|
523 |
#ifdef ALLOW_SIHEFF_CONTROL |
524 |
else if ( grdchkvarindex .eq. 42 ) then |
525 |
il=ilnblnk( xx_siheff_file ) |
526 |
write(fname(1:80),'(80a)') ' ' |
527 |
write(fname(1:80),'(3a,i10.10)') |
528 |
& yadmark, xx_siheff_file(1:il),'.',optimcycle |
529 |
call active_read_xy( fname, tmpfld2d, icvrec, |
530 |
& doglobalread, ladinit, optimcycle, |
531 |
& mythid, dummy) |
532 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
533 |
#endif /* ALLOW_SIHEFF_CONTROL */ |
534 |
|
535 |
#ifdef ALLOW_SIHSNOW_CONTROL |
536 |
else if ( grdchkvarindex .eq. 43 ) then |
537 |
il=ilnblnk( xx_sihsnow_file ) |
538 |
write(fname(1:80),'(80a)') ' ' |
539 |
write(fname(1:80),'(3a,i10.10)') |
540 |
& yadmark, xx_sihsnow_file(1:il),'.',optimcycle |
541 |
call active_read_xy( fname, tmpfld2d, icvrec, |
542 |
& doglobalread, ladinit, optimcycle, |
543 |
& mythid, dummy) |
544 |
xx_comp = tmpfld2d( itilepos,jtilepos,itile,jtile ) |
545 |
#endif /* ALLOW_SIHSNOW_CONTROL */ |
546 |
|
547 |
else |
548 |
ce --> this index does not exist yet. |
549 |
endif |
550 |
|
551 |
#endif /* ALLOW_GRDCHK */ |
552 |
|
553 |
end |
554 |
|