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

Annotation of /MITgcm/pkg/grdchk/grdchk_getxx.F

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


Revision 1.8 - (hide annotations) (download)
Tue Jun 24 16:08:45 2003 UTC (20 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51, checkpoint51b_pre, checkpoint51b_post, checkpoint51c_post, checkpoint51a_post
Changes since 1.7: +101 -1 lines
Merging for c51 vs. e34

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

  ViewVC Help
Powered by ViewVC 1.1.22