/[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.10 - (show annotations) (download)
Mon Oct 27 22:32:55 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint57g_pre, checkpoint57b_post, checkpoint52d_pre, checkpoint57g_post, checkpoint56b_post, checkpoint52j_pre, checkpoint54d_post, checkpoint54e_post, checkpoint57d_post, checkpoint57i_post, checkpoint52l_post, checkpoint52k_post, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint52f_post, checkpoint54f_post, checkpoint51t_post, checkpoint55i_post, checkpoint57l_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint55c_post, checkpoint52e_pre, checkpoint57f_post, checkpoint52e_post, checkpoint53d_post, checkpoint57a_post, checkpoint57h_pre, checkpoint52b_pre, checkpoint54b_post, checkpoint57h_post, checkpoint52m_post, checkpoint55g_post, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint57c_post, checkpoint52f_pre, checkpoint55d_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint57c_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint51r_post, checkpoint57e_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint52d_post, checkpoint53g_post, eckpoint57e_pre, checkpoint52a_pre, checkpoint52i_post, checkpoint52h_pre, checkpoint56a_post, checkpoint53f_post, checkpoint57h_done, checkpoint52j_post, checkpoint57j_post, checkpoint57f_pre, branch-netcdf, checkpoint52n_post, checkpoint53b_pre, checkpoint56c_post, checkpoint57a_pre, checkpoint55a_post, checkpoint51o_post, checkpoint57k_post, checkpoint53b_post, checkpoint52a_post, ecco_c52_e35, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, netcdf-sm0
Changes since 1.9: +3 -3 lines
o cleaning ALLOW_GRADIENT_CHECK -> ALLOW_GRDCHK
o cleaning some ALLOW_TANGENTLINEAR_RUN -> ALLOW_AUTODIFF
o bug fix in find_alpha.F for MDJWF:
  - modif. to alpha = 1/D*( dN/dT - rho*dD/Dt) to account for
    change rho -> rho-rhoConst
  - replace call find_rho to find_rhonum

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

  ViewVC Help
Powered by ViewVC 1.1.22