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

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

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


Revision 1.25 - (show annotations) (download)
Fri Oct 16 19:12:09 2009 UTC (14 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62, checkpoint62b, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.24: +2 -2 lines
Bug fix for new ctrl gen2d

1 C $Header: /u/gcmpack/MITgcm/pkg/grdchk/grdchk_setxx.F,v 1.24 2009/10/14 20:10:13 heimbach Exp $
2 C $Name: $
3
4 #include "CTRL_CPPOPTIONS.h"
5
6
7 subroutine grdchk_setxx(
8 I icvrec,
9 I theSimulationMode,
10 I itile,
11 I jtile,
12 I layer,
13 I itilepos,
14 I jtilepos,
15 I xx_comp_ref,
16 I mythid
17 & )
18
19 c ==================================================================
20 c SUBROUTINE grdchk_setxx
21 c ==================================================================
22 c
23 c o Set component a component of the control vector; xx(loc)
24 c
25 c started: Christian Eckert eckert@mit.edu 08-Mar-2000
26 c continued: heimbach@mit.edu: 13-Jun-2001
27 c
28 c ==================================================================
29 c SUBROUTINE grdchk_setxx
30 c ==================================================================
31
32 implicit none
33
34 c == global variables ==
35
36 #include "EEPARAMS.h"
37 #include "SIZE.h"
38 #include "ctrl.h"
39 #include "optim.h"
40 #include "grdchk.h"
41
42 c == routine arguments ==
43
44 integer icvrec
45 integer theSimulationMode
46 integer jtile
47 integer itile
48 integer layer
49 integer itilepos
50 integer jtilepos
51 _RL xx_comp_ref
52 integer mythid
53
54 #ifdef ALLOW_GRDCHK
55 c == local variables ==
56
57 integer i,j,k
58 integer il
59 integer dumiter
60 _RL dumtime
61 _RL dummy
62
63 logical doglobalread
64 logical ladinit
65
66 character*(80) fname
67
68 c-- == external ==
69
70 integer ilnblnk
71 external ilnblnk
72
73 c-- == end of interface ==
74
75 doglobalread = .false.
76 ladinit = .false.
77 dumiter = 0
78 dumtime = 0. _d 0
79 write(fname(1:80),'(80a)') ' '
80
81 if ( grdchkvarindex .eq. 0 ) then
82 STOP 'GRDCHK INDEX 0 NOT ALLOWED'
83
84 #ifdef ALLOW_THETA0_CONTROL
85 else if ( grdchkvarindex .eq. 1 ) then
86 il=ilnblnk( xx_theta_file )
87 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
88 write(fname(1:80),'(3a,i10.10)')
89 & yadmark, xx_theta_file(1:il),'.',optimcycle
90 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
91 write(fname(1:80),'(2a,i10.10)')
92 & xx_theta_file(1:il),'.',optimcycle
93 end if
94 call active_read_xyz( fname, tmpfld3d, 1,
95 & doglobalread, ladinit, optimcycle,
96 & mythid, dummy)
97 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
98 call active_write_xyz( fname, tmpfld3d, 1,
99 & optimcycle,
100 & mythid, dummy)
101 #endif /* ALLOW_THETA0_CONTROL */
102
103 #ifdef ALLOW_SALT0_CONTROL
104 else if ( grdchkvarindex .eq. 2 ) then
105 il=ilnblnk( xx_salt_file )
106 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
107 write(fname(1:80),'(3a,i10.10)')
108 & yadmark, xx_salt_file(1:il),'.',optimcycle
109 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
110 write(fname(1:80),'(2a,i10.10)')
111 & xx_salt_file(1:il),'.',optimcycle
112 end if
113 call active_read_xyz( fname, tmpfld3d, 1,
114 & doglobalread, ladinit, optimcycle,
115 & mythid, dummy)
116 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
117 call active_write_xyz( fname, tmpfld3d, 1,
118 & optimcycle,
119 & mythid, dummy)
120 #endif /* ALLOW_SALT0_CONTROL */
121
122 #ifdef ALLOW_HFLUX_CONTROL
123 else if ( grdchkvarindex .eq. 3 ) then
124 il=ilnblnk( xx_hflux_file )
125 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
126 write(fname(1:80),'(3a,i10.10)')
127 & yadmark, xx_hflux_file(1:il),'.',optimcycle
128 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
129 write(fname(1:80),'(2a,i10.10)')
130 & xx_hflux_file(1:il),'.',optimcycle
131 end if
132 call active_read_xy( fname, tmpfld2d, icvrec,
133 & doglobalread, ladinit, optimcycle,
134 & mythid, dummy)
135 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
136 call active_write_xy( fname, tmpfld2d, icvrec,
137 & optimcycle,
138 & mythid, dummy)
139 #endif /* ALLOW_HFLUX_CONTROL */
140
141 #ifdef ALLOW_SFLUX_CONTROL
142 else if ( grdchkvarindex .eq. 4 ) then
143 il=ilnblnk( xx_sflux_file )
144 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
145 write(fname(1:80),'(3a,i10.10)')
146 & yadmark, xx_sflux_file(1:il),'.',optimcycle
147 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
148 write(fname(1:80),'(2a,i10.10)')
149 & xx_sflux_file(1:il),'.',optimcycle
150 end if
151 call active_read_xy( fname, tmpfld2d, icvrec,
152 & doglobalread, ladinit, optimcycle,
153 & mythid, dummy)
154 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
155 call active_write_xy( fname, tmpfld2d, icvrec,
156 & optimcycle,
157 & mythid, dummy)
158 #endif /* ALLOW_SFLUX_CONTROL */
159
160 #ifdef ALLOW_USTRESS_CONTROL
161 else if ( grdchkvarindex .eq. 5 ) then
162 il=ilnblnk( xx_tauu_file )
163 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
164 write(fname(1:80),'(3a,i10.10)')
165 & yadmark, xx_tauu_file(1:il),'.',optimcycle
166 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
167 write(fname(1:80),'(2a,i10.10)')
168 & xx_tauu_file(1:il),'.',optimcycle
169 end if
170 call active_read_xy( fname, tmpfld2d, icvrec,
171 & doglobalread, ladinit, optimcycle,
172 & mythid, dummy)
173 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
174 call active_write_xy( fname, tmpfld2d, icvrec,
175 & optimcycle,
176 & mythid, dummy)
177 #endif /* ALLOW_USTRESS_CONTROL */
178
179 #ifdef ALLOW_VSTRESS_CONTROL
180 else if ( grdchkvarindex .eq. 6 ) then
181 il=ilnblnk( xx_tauv_file )
182 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
183 write(fname(1:80),'(3a,i10.10)')
184 & yadmark, xx_tauv_file(1:il),'.',optimcycle
185 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
186 write(fname(1:80),'(2a,i10.10)')
187 & xx_tauv_file(1:il),'.',optimcycle
188 end if
189 call active_read_xy( fname, tmpfld2d, icvrec,
190 & doglobalread, ladinit, optimcycle,
191 & mythid, dummy)
192 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
193 call active_write_xy( fname, tmpfld2d, icvrec,
194 & optimcycle,
195 & mythid, dummy)
196 #endif /* ALLOW_VSTRESS_CONTROL */
197
198 #ifdef ALLOW_ATEMP_CONTROL
199 else if ( grdchkvarindex .eq. 7 ) then
200 il=ilnblnk( xx_atemp_file )
201 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
202 write(fname(1:80),'(3a,i10.10)')
203 & yadmark, xx_atemp_file(1:il),'.',optimcycle
204 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
205 write(fname(1:80),'(2a,i10.10)')
206 & xx_atemp_file(1:il),'.',optimcycle
207 end if
208 call active_read_xy( fname, tmpfld2d, icvrec,
209 & doglobalread, ladinit, optimcycle,
210 & mythid, dummy)
211 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
212 call active_write_xy( fname, tmpfld2d, icvrec,
213 & optimcycle,
214 & mythid, dummy)
215 #endif /* ALLOW_ATEMP_CONTROL */
216
217 #ifdef ALLOW_AQH_CONTROL
218 else if ( grdchkvarindex .eq. 8 ) then
219 il=ilnblnk( xx_aqh_file )
220 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
221 write(fname(1:80),'(3a,i10.10)')
222 & yadmark, xx_aqh_file(1:il),'.',optimcycle
223 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
224 write(fname(1:80),'(2a,i10.10)')
225 & xx_aqh_file(1:il),'.',optimcycle
226 end if
227 call active_read_xy( fname, tmpfld2d, icvrec,
228 & doglobalread, ladinit, optimcycle,
229 & mythid, dummy)
230 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
231 call active_write_xy( fname, tmpfld2d, icvrec,
232 & optimcycle,
233 & mythid, dummy)
234 #endif /* ALLOW_AQH_CONTROL */
235
236 #ifdef ALLOW_UWIND_CONTROL
237 else if ( grdchkvarindex .eq. 9 ) then
238 il=ilnblnk( xx_uwind_file )
239 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
240 write(fname(1:80),'(3a,i10.10)')
241 & yadmark, xx_uwind_file(1:il),'.',optimcycle
242 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
243 write(fname(1:80),'(2a,i10.10)')
244 & xx_uwind_file(1:il),'.',optimcycle
245 end if
246 call active_read_xy( fname, tmpfld2d, icvrec,
247 & doglobalread, ladinit, optimcycle,
248 & mythid, dummy)
249 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
250 call active_write_xy( fname, tmpfld2d, icvrec,
251 & optimcycle,
252 & mythid, dummy)
253 #endif /* ALLOW_UWIND_CONTROL */
254
255 #ifdef ALLOW_VWIND_CONTROL
256 else if ( grdchkvarindex .eq. 10 ) then
257 il=ilnblnk( xx_vwind_file )
258 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
259 write(fname(1:80),'(3a,i10.10)')
260 & yadmark, xx_vwind_file(1:il),'.',optimcycle
261 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
262 write(fname(1:80),'(2a,i10.10)')
263 & xx_vwind_file(1:il),'.',optimcycle
264 end if
265 call active_read_xy( fname, tmpfld2d, icvrec,
266 & doglobalread, ladinit, optimcycle,
267 & mythid, dummy)
268 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
269 call active_write_xy( fname, tmpfld2d, icvrec,
270 & optimcycle,
271 & mythid, dummy)
272 #endif /* ALLOW_VWIND_CONTROL */
273
274 #ifdef ALLOW_OBCSN_CONTROL
275 else if ( grdchkvarindex .eq. 11 ) then
276 il=ilnblnk( xx_obcsn_file )
277 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
278 write(fname(1:80),'(3a,i10.10)')
279 & yadmark, xx_obcsn_file(1:il),'.',optimcycle
280 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
281 write(fname(1:80),'(2a,i10.10)')
282 & xx_obcsn_file(1:il),'.',optimcycle
283 end if
284
285 call active_read_xz( fname, tmpfldxz, icvrec,
286 & doglobalread, ladinit, optimcycle,
287 & mythid, dummy)
288
289 tmpfldxz( itilepos,layer,itile,jtile ) = xx_comp_ref
290
291 call active_write_xz( fname, tmpfldxz, icvrec,
292 & optimcycle,
293 & mythid, dummy)
294
295 #endif /* ALLOW_OBCSN_CONTROL */
296
297 #ifdef ALLOW_OBCSS_CONTROL
298 else if ( grdchkvarindex .eq. 12 ) then
299 il=ilnblnk( xx_obcss_file )
300 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
301 write(fname(1:80),'(3a,i10.10)')
302 & yadmark, xx_obcss_file(1:il),'.',optimcycle
303 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
304 write(fname(1:80),'(2a,i10.10)')
305 & xx_obcss_file(1:il),'.',optimcycle
306 end if
307
308 call active_read_xz( fname, tmpfldxz, icvrec,
309 & doglobalread, ladinit, optimcycle,
310 & mythid, dummy)
311
312 tmpfldxz( itilepos,layer,itile,jtile ) = xx_comp_ref
313
314 call active_write_xz( fname, tmpfldxz, icvrec,
315 & optimcycle,
316 & mythid, dummy)
317
318 #endif /* ALLOW_OBCSS_CONTROL */
319
320 #ifdef ALLOW_OBCSW_CONTROL
321 else if ( grdchkvarindex .eq. 13 ) then
322 il=ilnblnk( xx_obcsw_file )
323 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
324 write(fname(1:80),'(3a,i10.10)')
325 & yadmark, xx_obcsw_file(1:il),'.',optimcycle
326 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
327 write(fname(1:80),'(2a,i10.10)')
328 & xx_obcsw_file(1:il),'.',optimcycle
329 end if
330
331 call active_read_yz( fname, tmpfldyz, icvrec,
332 & doglobalread, ladinit, optimcycle,
333 & mythid, dummy)
334
335 tmpfldyz( jtilepos,layer,itile,jtile ) = xx_comp_ref
336
337 call active_write_yz( fname, tmpfldyz, icvrec,
338 & optimcycle,
339 & mythid, dummy)
340
341 #endif /* ALLOW_OBCSW_CONTROL */
342
343 #ifdef ALLOW_OBCSE_CONTROL
344 else if ( grdchkvarindex .eq. 14 ) then
345 il=ilnblnk( xx_obcse_file )
346 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
347 write(fname(1:80),'(3a,i10.10)')
348 & yadmark, xx_obcse_file(1:il),'.',optimcycle
349 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
350 write(fname(1:80),'(2a,i10.10)')
351 & xx_obcse_file(1:il),'.',optimcycle
352 end if
353
354 call active_read_yz( fname, tmpfldyz, icvrec,
355 & doglobalread, ladinit, optimcycle,
356 & mythid, dummy)
357
358 tmpfldyz( jtilepos,layer,itile,jtile ) = xx_comp_ref
359
360 call active_write_yz( fname, tmpfldyz, icvrec,
361 & optimcycle,
362 & mythid, dummy)
363
364 #endif /* ALLOW_OBCSE_CONTROL */
365
366 #ifdef ALLOW_DIFFKR_CONTROL
367 else if ( grdchkvarindex .eq. 15 ) then
368 il=ilnblnk( xx_diffkr_file )
369 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
370 write(fname(1:80),'(3a,i10.10)')
371 & yadmark, xx_diffkr_file(1:il),'.',optimcycle
372 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
373 write(fname(1:80),'(2a,i10.10)')
374 & xx_diffkr_file(1:il),'.',optimcycle
375 end if
376
377 call active_read_xyz( fname, tmpfld3d, 1,
378 & doglobalread, ladinit, optimcycle,
379 & mythid, dummy)
380
381 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
382
383 call active_write_xyz( fname, tmpfld3d, 1,
384 & optimcycle,
385 & mythid, dummy)
386
387 #endif /* ALLOW_DIFFKR_CONTROL */
388
389 #ifdef ALLOW_KAPGM_CONTROL
390 else if ( grdchkvarindex .eq. 16 ) then
391 il=ilnblnk( xx_kapgm_file )
392 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
393 write(fname(1:80),'(3a,i10.10)')
394 & yadmark, xx_kapgm_file(1:il),'.',optimcycle
395 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
396 write(fname(1:80),'(2a,i10.10)')
397 & xx_kapgm_file(1:il),'.',optimcycle
398 end if
399
400 call active_read_xyz( fname, tmpfld3d, 1,
401 & doglobalread, ladinit, optimcycle,
402 & mythid, dummy)
403
404 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
405
406 call active_write_xyz( fname, tmpfld3d, 1,
407 & optimcycle,
408 & mythid, dummy)
409
410 #endif /* ALLOW_KAPGM_CONTROL */
411
412 #ifdef ALLOW_KAPREDI_CONTROL
413 else if ( grdchkvarindex .eq. 16 ) then
414 il=ilnblnk( xx_kapredi_file )
415 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
416 write(fname(1:80),'(3a,i10.10)')
417 & yadmark, xx_kapredi_file(1:il),'.',optimcycle
418 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
419 write(fname(1:80),'(2a,i10.10)')
420 & xx_kapredi_file(1:il),'.',optimcycle
421 end if
422
423 call active_read_xyz( fname, tmpfld3d, 1,
424 & doglobalread, ladinit, optimcycle,
425 & mythid, dummy)
426
427 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
428
429 call active_write_xyz( fname, tmpfld3d, 1,
430 & optimcycle,
431 & mythid, dummy)
432
433 #endif /* ALLOW_KAPREDI_CONTROL */
434
435 #ifdef ALLOW_TR10_CONTROL
436 else if ( grdchkvarindex .eq. 17 ) then
437 il=ilnblnk( xx_tr1_file )
438 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
439 write(fname(1:80),'(3a,i10.10)')
440 & yadmark, xx_tr1_file(1:il),'.',optimcycle
441 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
442 write(fname(1:80),'(2a,i10.10)')
443 & xx_tr1_file(1:il),'.',optimcycle
444 end if
445
446 call active_read_xyz( fname, tmpfld3d, 1,
447 & doglobalread, ladinit, optimcycle,
448 & mythid, dummy)
449
450 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
451
452 call active_write_xyz( fname, tmpfld3d, 1,
453 & optimcycle,
454 & mythid, dummy)
455
456 #endif /* ALLOW_TR10_CONTROL */
457
458 #if (defined (ALLOW_SST_CONTROL) || defined (ALLOW_SST0_CONTROL))
459 else if ( grdchkvarindex .eq. 18 ) then
460 il=ilnblnk( xx_sst_file )
461 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
462 write(fname(1:80),'(3a,i10.10)')
463 & yadmark, xx_sst_file(1:il),'.',optimcycle
464 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
465 write(fname(1:80),'(2a,i10.10)')
466 & xx_sst_file(1:il),'.',optimcycle
467 end if
468 call active_read_xy( fname, tmpfld2d, icvrec,
469 & doglobalread, ladinit, optimcycle,
470 & mythid, dummy)
471 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
472 call active_write_xy( fname, tmpfld2d, icvrec,
473 & optimcycle,
474 & mythid, dummy)
475 #endif /* ALLOW_SST0_CONTROL */
476
477 #if (defined (ALLOW_SSS_CONTROL) || defined (ALLOW_SSS0_CONTROL))
478 else if ( grdchkvarindex .eq. 19 ) then
479 il=ilnblnk( xx_sss_file )
480 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
481 write(fname(1:80),'(3a,i10.10)')
482 & yadmark, xx_sss_file(1:il),'.',optimcycle
483 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
484 write(fname(1:80),'(2a,i10.10)')
485 & xx_sss_file(1:il),'.',optimcycle
486 end if
487 call active_read_xy( fname, tmpfld2d, icvrec,
488 & doglobalread, ladinit, optimcycle,
489 & mythid, dummy)
490 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
491 call active_write_xy( fname, tmpfld2d, icvrec,
492 & optimcycle,
493 & mythid, dummy)
494 #endif /* ALLOW_SSS0_CONTROL */
495
496 #ifdef ALLOW_DEPTH_CONTROL
497 else if ( grdchkvarindex .eq. 20 ) then
498 il=ilnblnk( xx_depth_file )
499 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
500 write(fname(1:80),'(3a,i10.10)')
501 & yadmark, xx_depth_file(1:il),'.',optimcycle
502 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
503 write(fname(1:80),'(2a,i10.10)')
504 & xx_depth_file(1:il),'.',optimcycle
505 end if
506 call active_read_xy( fname, tmpfld2d, icvrec,
507 & doglobalread, ladinit, optimcycle,
508 & mythid, dummy)
509 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
510 call active_write_xy( fname, tmpfld2d, icvrec,
511 & optimcycle,
512 & mythid, dummy)
513 #endif /* ALLOW_DEPTH_CONTROL */
514
515 #ifdef ALLOW_EFLUXY0_CONTROL
516 else if ( grdchkvarindex .eq. 21 ) then
517 il=ilnblnk( xx_efluxy_file )
518 write(fname(1:80),'(80a)') ' '
519 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
520 write(fname(1:80),'(3a,i10.10)')
521 & yadmark, xx_efluxy_file(1:il),'.',optimcycle
522 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
523 write(fname(1:80),'(2a,i10.10)')
524 & xx_efluxy_file(1:il),'.',optimcycle
525 end if
526 call active_read_xyz( fname, tmpfld3d, 1,
527 & doglobalread, ladinit, optimcycle,
528 & mythid, dummy)
529 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
530 call active_write_xyz( fname, tmpfld3d, 1,
531 & optimcycle,
532 & mythid, dummy)
533 #endif /* ALLOW_EFLUXY0_CONTROL */
534
535 #ifdef ALLOW_EFLUXP0_CONTROL
536 else if ( grdchkvarindex .eq. 22 ) then
537 il=ilnblnk( xx_efluxp_file )
538 write(fname(1:80),'(80a)') ' '
539 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
540 write(fname(1:80),'(3a,i10.10)')
541 & yadmark, xx_efluxp_file(1:il),'.',optimcycle
542 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
543 write(fname(1:80),'(2a,i10.10)')
544 & xx_efluxp_file(1:il),'.',optimcycle
545 end if
546 call active_read_xyz( fname, tmpfld3d, 1,
547 & doglobalread, ladinit, optimcycle,
548 & mythid, dummy)
549 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
550 call active_write_xyz( fname, tmpfld3d, 1,
551 & optimcycle,
552 & mythid, dummy)
553 #endif /* ALLOW_EFLUXP0_CONTROL */
554
555 #ifdef ALLOW_HFLUXM_CONTROL
556 else if ( grdchkvarindex .eq. 24 ) then
557 il=ilnblnk( xx_hfluxm_file )
558 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
559 write(fname(1:80),'(3a,i10.10)')
560 & yadmark, xx_hfluxm_file(1:il),'.',optimcycle
561 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
562 write(fname(1:80),'(2a,i10.10)')
563 & xx_hfluxm_file(1:il),'.',optimcycle
564 end if
565 call active_read_xy( fname, tmpfld2d, icvrec,
566 & doglobalread, ladinit, optimcycle,
567 & mythid, dummy)
568 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
569 call active_write_xy( fname, tmpfld2d, icvrec,
570 & optimcycle,
571 & mythid, dummy)
572 #endif /* ALLOW_HFLUXM_CONTROL */
573
574 #ifdef ALLOW_GEN2D_CONTROL
575 else if ( grdchkvarindex .eq. 30 ) then
576 il=ilnblnk( xx_gen2d_file )
577 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
578 write(fname(1:80),'(3a,i10.10)')
579 & yadmark, xx_gen2d_file(1:il),'.',optimcycle
580 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
581 write(fname(1:80),'(2a,i10.10)')
582 & xx_gen2d_file(1:il),'.',optimcycle
583 end if
584 call active_read_xy( fname, tmpfld2d, icvrec,
585 & doglobalread, ladinit, optimcycle,
586 & mythid, dummy)
587 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
588 call active_write_xy( fname, tmpfld2d, icvrec,
589 & optimcycle,
590 & mythid, dummy)
591 #endif /* ALLOW_GEN2D_CONTROL */
592
593 #ifdef ALLOW_GEN3D_CONTROL
594 else if ( grdchkvarindex .eq. 31 ) then
595 il=ilnblnk( xx_gen3d_file )
596 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
597 write(fname(1:80),'(3a,i10.10)')
598 & yadmark, xx_gen3d_file(1:il),'.',optimcycle
599 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
600 write(fname(1:80),'(2a,i10.10)')
601 & xx_gen3d_file(1:il),'.',optimcycle
602 end if
603
604 call active_read_xyz( fname, tmpfld3d, 1,
605 & doglobalread, ladinit, optimcycle,
606 & mythid, dummy)
607
608 tmpfld3d( itilepos,jtilepos,layer,itile,jtile ) = xx_comp_ref
609
610 call active_write_xyz( fname, tmpfld3d, 1,
611 & optimcycle,
612 & mythid, dummy)
613 #endif /* ALLOW_GEN3D_CONTROL */
614
615 #ifdef ALLOW_PRECIP_CONTROL
616 else if ( grdchkvarindex .eq. 32 ) then
617 il=ilnblnk( xx_precip_file )
618 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
619 write(fname(1:80),'(3a,i10.10)')
620 & yadmark, xx_precip_file(1:il),'.',optimcycle
621 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
622 write(fname(1:80),'(2a,i10.10)')
623 & xx_precip_file(1:il),'.',optimcycle
624 end if
625 call active_read_xy( fname, tmpfld2d, icvrec,
626 & doglobalread, ladinit, optimcycle,
627 & mythid, dummy)
628 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
629 call active_write_xy( fname, tmpfld2d, icvrec,
630 & optimcycle,
631 & mythid, dummy)
632 #endif /* ALLOW_PRECIP_CONTROL */
633
634 #ifdef ALLOW_SWFLUX_CONTROL
635 else if ( grdchkvarindex .eq. 33 ) then
636 il=ilnblnk( xx_swflux_file )
637 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
638 write(fname(1:80),'(3a,i10.10)')
639 & yadmark, xx_swflux_file(1:il),'.',optimcycle
640 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
641 write(fname(1:80),'(2a,i10.10)')
642 & xx_swflux_file(1:il),'.',optimcycle
643 end if
644 call active_read_xy( fname, tmpfld2d, icvrec,
645 & doglobalread, ladinit, optimcycle,
646 & mythid, dummy)
647 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
648 call active_write_xy( fname, tmpfld2d, icvrec,
649 & optimcycle,
650 & mythid, dummy)
651 #endif /* ALLOW_SWFLUX_CONTROL */
652
653 #ifdef ALLOW_SWDOWN_CONTROL
654 else if ( grdchkvarindex .eq. 34 ) then
655 il=ilnblnk( xx_swdown_file )
656 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
657 write(fname(1:80),'(3a,i10.10)')
658 & yadmark, xx_swdown_file(1:il),'.',optimcycle
659 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
660 write(fname(1:80),'(2a,i10.10)')
661 & xx_swdown_file(1:il),'.',optimcycle
662 end if
663 call active_read_xy( fname, tmpfld2d, icvrec,
664 & doglobalread, ladinit, optimcycle,
665 & mythid, dummy)
666 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
667 call active_write_xy( fname, tmpfld2d, icvrec,
668 & optimcycle,
669 & mythid, dummy)
670 #endif /* ALLOW_SWDOWN_CONTROL */
671
672 #ifdef ALLOW_LWFLUX_CONTROL
673 else if ( grdchkvarindex .eq. 35 ) then
674 il=ilnblnk( xx_lwflux_file )
675 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
676 write(fname(1:80),'(3a,i10.10)')
677 & yadmark, xx_lwflux_file(1:il),'.',optimcycle
678 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
679 write(fname(1:80),'(2a,i10.10)')
680 & xx_lwflux_file(1:il),'.',optimcycle
681 end if
682 call active_read_xy( fname, tmpfld2d, icvrec,
683 & doglobalread, ladinit, optimcycle,
684 & mythid, dummy)
685 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
686 call active_write_xy( fname, tmpfld2d, icvrec,
687 & optimcycle,
688 & mythid, dummy)
689 #endif /* ALLOW_LWFLUX_CONTROL */
690
691 #ifdef ALLOW_LWDOWN_CONTROL
692 else if ( grdchkvarindex .eq. 36 ) then
693 il=ilnblnk( xx_lwdown_file )
694 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
695 write(fname(1:80),'(3a,i10.10)')
696 & yadmark, xx_lwdown_file(1:il),'.',optimcycle
697 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
698 write(fname(1:80),'(2a,i10.10)')
699 & xx_lwdown_file(1:il),'.',optimcycle
700 end if
701 call active_read_xy( fname, tmpfld2d, icvrec,
702 & doglobalread, ladinit, optimcycle,
703 & mythid, dummy)
704 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
705 call active_write_xy( fname, tmpfld2d, icvrec,
706 & optimcycle,
707 & mythid, dummy)
708 #endif /* ALLOW_LWDOWN_CONTROL */
709
710 #ifdef ALLOW_EVAP_CONTROL
711 else if ( grdchkvarindex .eq. 37 ) then
712 il=ilnblnk( xx_evap_file )
713 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
714 write(fname(1:80),'(3a,i10.10)')
715 & yadmark, xx_evap_file(1:il),'.',optimcycle
716 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
717 write(fname(1:80),'(2a,i10.10)')
718 & xx_evap_file(1:il),'.',optimcycle
719 end if
720 call active_read_xy( fname, tmpfld2d, icvrec,
721 & doglobalread, ladinit, optimcycle,
722 & mythid, dummy)
723 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
724 call active_write_xy( fname, tmpfld2d, icvrec,
725 & optimcycle,
726 & mythid, dummy)
727 #endif /* ALLOW_EVAP_CONTROL */
728
729 #ifdef ALLOW_SNOWPRECIP_CONTROL
730 else if ( grdchkvarindex .eq. 38 ) then
731 il=ilnblnk( xx_snowprecip_file )
732 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
733 write(fname(1:80),'(3a,i10.10)')
734 & yadmark, xx_snowprecip_file(1:il),'.',optimcycle
735 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
736 write(fname(1:80),'(2a,i10.10)')
737 & xx_snowprecip_file(1:il),'.',optimcycle
738 end if
739 call active_read_xy( fname, tmpfld2d, icvrec,
740 & doglobalread, ladinit, optimcycle,
741 & mythid, dummy)
742 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
743 call active_write_xy( fname, tmpfld2d, icvrec,
744 & optimcycle,
745 & mythid, dummy)
746 #endif /* ALLOW_SNOWPRECIP_CONTROL */
747
748 #ifdef ALLOW_APRESSURE_CONTROL
749 else if ( grdchkvarindex .eq. 39 ) then
750 il=ilnblnk( xx_apressure_file )
751 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
752 write(fname(1:80),'(3a,i10.10)')
753 & yadmark, xx_apressure_file(1:il),'.',optimcycle
754 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
755 write(fname(1:80),'(2a,i10.10)')
756 & xx_apressure_file(1:il),'.',optimcycle
757 end if
758 call active_read_xy( fname, tmpfld2d, icvrec,
759 & doglobalread, ladinit, optimcycle,
760 & mythid, dummy)
761 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
762 call active_write_xy( fname, tmpfld2d, icvrec,
763 & optimcycle,
764 & mythid, dummy)
765 #endif /* ALLOW_APRESSURE_CONTROL */
766
767 #ifdef ALLOW_RUNOFF_CONTROL
768 else if ( grdchkvarindex .eq. 40 ) then
769 il=ilnblnk( xx_runoff_file )
770 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
771 write(fname(1:80),'(3a,i10.10)')
772 & yadmark, xx_runoff_file(1:il),'.',optimcycle
773 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
774 write(fname(1:80),'(2a,i10.10)')
775 & xx_runoff_file(1:il),'.',optimcycle
776 end if
777 call active_read_xy( fname, tmpfld2d, icvrec,
778 & doglobalread, ladinit, optimcycle,
779 & mythid, dummy)
780 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
781 call active_write_xy( fname, tmpfld2d, icvrec,
782 & optimcycle,
783 & mythid, dummy)
784 #endif /* ALLOW_RUNOFF_CONTROL */
785
786 #ifdef ALLOW_SIAREA_CONTROL
787 else if ( grdchkvarindex .eq. 41 ) then
788 il=ilnblnk( xx_siarea_file )
789 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
790 write(fname(1:80),'(3a,i10.10)')
791 & yadmark, xx_siarea_file(1:il),'.',optimcycle
792 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
793 write(fname(1:80),'(2a,i10.10)')
794 & xx_siarea_file(1:il),'.',optimcycle
795 end if
796 call active_read_xy( fname, tmpfld2d, icvrec,
797 & doglobalread, ladinit, optimcycle,
798 & mythid, dummy)
799 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
800 call active_write_xy( fname, tmpfld2d, icvrec,
801 & optimcycle,
802 & mythid, dummy)
803 #endif /* ALLOW_SIAREA_CONTROL */
804
805 #ifdef ALLOW_SIHEFF_CONTROL
806 else if ( grdchkvarindex .eq. 42 ) then
807 il=ilnblnk( xx_siheff_file )
808 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
809 write(fname(1:80),'(3a,i10.10)')
810 & yadmark, xx_siheff_file(1:il),'.',optimcycle
811 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
812 write(fname(1:80),'(2a,i10.10)')
813 & xx_siheff_file(1:il),'.',optimcycle
814 end if
815 call active_read_xy( fname, tmpfld2d, icvrec,
816 & doglobalread, ladinit, optimcycle,
817 & mythid, dummy)
818 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
819 call active_write_xy( fname, tmpfld2d, icvrec,
820 & optimcycle,
821 & mythid, dummy)
822 #endif /* ALLOW_SIHEFF_CONTROL */
823
824 #ifdef ALLOW_SIHSNOW_CONTROL
825 else if ( grdchkvarindex .eq. 43 ) then
826 il=ilnblnk( xx_sihsnow_file )
827 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
828 write(fname(1:80),'(3a,i10.10)')
829 & yadmark, xx_sihsnow_file(1:il),'.',optimcycle
830 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
831 write(fname(1:80),'(2a,i10.10)')
832 & xx_sihsnow_file(1:il),'.',optimcycle
833 end if
834 call active_read_xy( fname, tmpfld2d, icvrec,
835 & doglobalread, ladinit, optimcycle,
836 & mythid, dummy)
837 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
838 call active_write_xy( fname, tmpfld2d, icvrec,
839 & optimcycle,
840 & mythid, dummy)
841 #endif /* ALLOW_SIHSNOW_CONTROL */
842
843
844 #ifdef ALLOW_ETAN0_CONTROL
845 else if ( grdchkvarindex .eq. 29 ) then
846 il=ilnblnk( xx_etan_file )
847 if ( theSimulationMode .EQ. TANGENT_SIMULATION ) then
848 write(fname(1:80),'(3a,i10.10)')
849 & yadmark, xx_etan_file(1:il),'.',optimcycle
850 else if ( theSimulationMode .EQ. FORWARD_SIMULATION ) then
851 write(fname(1:80),'(2a,i10.10)')
852 & xx_etan_file(1:il),'.',optimcycle
853 end if
854 call active_read_xy( fname, tmpfld2d, 1,
855 & doglobalread, ladinit, optimcycle,
856 & mythid, dummy)
857 tmpfld2d( itilepos,jtilepos,itile,jtile ) = xx_comp_ref
858 print*,"fname etc.",fname,itilepos,jtilepos,itile,jtile,
859 & tmpfld2d( itilepos,jtilepos,itile,jtile )
860 call active_write_xy( fname, tmpfld2d, 1,
861 & optimcycle,
862 & mythid, dummy)
863 #endif /* ALLOW_ETAN0_CONTROL */
864
865 else
866 ce --> this index does not exist yet.
867 endif
868
869 #endif /* ALLOW_GRDCHK */
870
871 end
872

  ViewVC Help
Powered by ViewVC 1.1.22