/[MITgcm]/MITgcm/pkg/autodiff/active_file_control.F
ViewVC logotype

Annotation of /MITgcm/pkg/autodiff/active_file_control.F

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


Revision 1.4 - (hide annotations) (download)
Mon Sep 16 18:11:58 2002 UTC (21 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint47e_post, checkpoint46l_post, checkpoint46g_pre, checkpoint47c_post, checkpoint50c_post, checkpoint46f_post, checkpoint48e_post, checkpoint50c_pre, checkpoint48i_post, checkpoint46l_pre, checkpoint51, checkpoint50, checkpoint50d_post, checkpoint50b_pre, checkpoint51f_post, checkpoint48b_post, checkpoint51d_post, checkpoint48c_pre, checkpoint47d_pre, c49_autodiff, checkpoint47a_post, checkpoint48d_pre, checkpoint51j_post, checkpoint47i_post, checkpoint47d_post, checkpoint48d_post, checkpoint48f_post, checkpoint46j_pre, checkpoint48h_post, checkpoint51b_pre, checkpoint47g_post, checkpoint46j_post, checkpoint51h_pre, checkpoint46k_post, checkpoint48a_post, checkpoint50f_post, checkpoint50a_post, checkpoint50f_pre, checkpoint47j_post, branch-exfmods-tag, branchpoint-genmake2, checkpoint48c_post, checkpoint51b_post, checkpoint51c_post, checkpoint47b_post, checkpoint46h_pre, checkpoint46m_post, checkpoint50g_post, checkpoint46g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51i_pre, checkpoint47f_post, checkpoint50e_post, checkpoint46i_post, checkpoint50d_pre, checkpoint51e_post, checkpoint47, checkpoint48, checkpoint49, checkpoint46h_post, checkpoint51f_pre, checkpoint48g_post, checkpoint47h_post, checkpoint51g_post, checkpoint50b_post, checkpoint51a_post
Branch point for: branch-exfmods-curt, branch-genmake2
Changes since 1.3: +130 -0 lines
Enable tangent linear (forward mode) gradient checks:
o extended active file handling to g_... files
o added TANGENT_SIMULATION to theSimulationMode
o extended grdchk package accordingly

1 heimbach 1.3
2     #include "CPP_OPTIONS.h"
3    
4     c ==================================================================
5     c
6     c active_file_control.F: Routines to handle the i/o of active vari-
7     c ables for the adjoint calculations. All
8     c files are direct access files.
9     c
10     c Routines:
11     c
12     c o active_read_rl - Basic routine to handle active read
13     c operations.
14     c o active_write_rl - Basic routine to handle active write
15     c operations.
16     c
17     c o active_read_tile_rl - Basic routine to handle active read
18     c operations.
19     c o active_write_tile_rl - Basic routine to handle active write
20     c operations.
21     c
22     c
23     c changed: Christian Eckert eckert@mit.edu 24-Apr-2000
24     c - Added routines that do active writes on tiles
25     c instead of a whole thread.
26     c
27     c ==================================================================
28    
29    
30     CBOP
31     C !ROUTINE: active_read_rl
32     C !INTERFACE:
33     subroutine active_read_rl(
34     I active_var_file,
35     O active_var,
36     I globalfile,
37     I lAdInit,
38     I irec,
39     I mynr,
40     I theSimulationMode,
41     I myOptimIter,
42     I mythid
43     & )
44    
45     C !DESCRIPTION: \bv
46     c ==================================================================
47     c SUBROUTINE active_read_rl
48     c ==================================================================
49     c o Read an active _RL variable from file.
50     c The variable *globalfile* can be used as a switch, which allows
51     c to read from a global file. The adjoint files are, however, always
52     c treated as tiled files.
53     c started: Christian Eckert eckert@mit.edu Jan-1999
54     c ==================================================================
55     c SUBROUTINE active_read_rl
56     c ==================================================================
57     C \ev
58    
59     C !USES:
60     implicit none
61    
62     c == global variables ==
63     #include "EEPARAMS.h"
64     #include "SIZE.h"
65     #include "PARAMS.h"
66    
67     C !INPUT/OUTPUT PARAMETERS:
68     c == routine arguments ==
69     c active_var_file: filename
70     c active_var: array
71     c irec: record number
72     c myOptimIter: number of optimization iteration (default: 0)
73     c mythid: thread number for this instance
74     c doglobalread: flag for global or local read/write
75     c (default: .false.)
76     c lAdInit: initialisation of corresponding adjoint
77     c variable and write to active file
78     c mynr: vertical array dimension
79     c theSimulationMode: forward mode or reverse mode simulation
80     character*(*) active_var_file
81     logical globalfile
82     logical lAdInit
83     integer irec
84     integer mynr
85     integer theSimulationMode
86     integer myOptimIter
87     integer mythid
88     _RL active_var(1-olx:snx+olx,1-oly:sny+oly,mynr,nsx,nsy)
89    
90     C !LOCAL VARIABLES:
91     c == local variables ==
92     character*(2) adpref
93     character*(80) adfname
94     integer bi,bj
95     integer i,j,k
96     integer oldprec
97     integer prec
98     integer il
99     integer ilnblnk
100     logical writeglobalfile
101     _RL active_data_t(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
102    
103     c == functions ==
104     external ilnblnk
105    
106     c == end of interface ==
107     CEOP
108    
109     c force 64-bit io
110     oldPrec = readBinaryPrec
111     readBinaryPrec = precFloat64
112     prec = precFloat64
113    
114     write(adfname(1:80),'(80a)') ' '
115     adpref = 'ad'
116     il = ilnblnk( active_var_file )
117    
118     write(adfname(1:2),'(a)') adpref
119     write(adfname(3:il+2),'(a)') active_var_file(1:il)
120    
121     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
122     c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<<
123     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
124    
125     if (theSimulationMode .eq. FORWARD_SIMULATION) then
126    
127     _BEGIN_MASTER( mythid )
128    
129     c Read the active variable from file.
130    
131     call mdsreadfield(
132     & active_var_file,
133     & prec,
134     & 'RL',
135     & mynr,
136     & active_var,
137     & irec,
138     & mythid )
139    
140     if (lAdInit) then
141     c Initialise the corresponding adjoint variable on the
142     c adjoint variable's file. These files are tiled.
143    
144     writeglobalfile = .false.
145     do bj = 1,nsy
146     do bi = 1,nsx
147     do j=1,sny
148     do i=1,snx
149     active_data_t(i,j,bi,bj)= 0. _d 0
150     enddo
151     enddo
152     enddo
153     enddo
154    
155     do k = 1,mynr
156     call mdswritefield(
157     & adfname,
158     & prec,
159     & globalfile,
160     & 'RL',
161     & 1,
162     & active_data_t,
163     & (irec-1)*mynr+k,
164     & myOptimIter,
165     & mythid )
166     enddo
167     endif
168    
169     _END_MASTER( mythid )
170    
171     endif
172    
173     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
174     c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<<
175     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
176    
177     if (theSimulationMode .eq. REVERSE_SIMULATION) then
178    
179     _BEGIN_MASTER( mythid )
180    
181     writeglobalfile = .false.
182     do k=1,mynr
183     c Read data from file layer by layer.
184     call mdsreadfield(
185     & active_var_file,
186     & prec,
187     & 'RL',
188     & 1,
189     & active_data_t,
190     & (irec-1)*mynr+k,
191     & mythid )
192    
193     c Add active_var from appropriate location to data.
194     do bj = 1,nsy
195     do bi = 1,nsx
196     do j=1,sny
197     do i=1,snx
198     active_data_t(i,j,bi,bj) = active_data_t(i,j,bi,bj) +
199     & active_var(i,j,k,bi,bj)
200     enddo
201     enddo
202     enddo
203     enddo
204    
205     c Store the result on disk.
206     call mdswritefield(
207     & active_var_file,
208     & prec,
209     & writeglobalfile,
210     & 'RL',
211     & 1,
212     & active_data_t,
213     & (irec-1)*mynr+k,
214     & myOptimIter,
215     & mythid )
216     enddo
217    
218    
219     c Set active_var to zero.
220     do k=1,mynr
221     do bj = 1,nsy
222     do bi = 1,nsx
223     do j=1,sny
224     do i=1,snx
225     active_var(i,j,k,bi,bj) = 0. _d 0
226     enddo
227     enddo
228     enddo
229     enddo
230     enddo
231    
232     _END_MASTER( mythid )
233     endif
234    
235 heimbach 1.4 c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
236     c >>>>>>>>>>>>>>>>>>> TANGENT RUN <<<<<<<<<<<<<<<<<<<
237     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
238    
239     if (theSimulationMode .eq. TANGENT_SIMULATION) then
240    
241     _BEGIN_MASTER( mythid )
242    
243     c Read the active variable from file.
244    
245     call mdsreadfield(
246     & active_var_file,
247     & prec,
248     & 'RL',
249     & mynr,
250     & active_var,
251     & irec,
252     & mythid )
253    
254     _END_MASTER( mythid )
255     endif
256    
257 heimbach 1.3 c Reset default io precision.
258     readBinaryPrec = oldPrec
259    
260     _BARRIER
261    
262     return
263     end
264    
265    
266     CBOP
267     C !ROUTINE: active_write_rl
268     C !INTERFACE:
269     subroutine active_write_rl(
270     I active_var_file,
271     I active_var,
272     I globalfile,
273     I irec,
274     I mynr,
275     I theSimulationMode,
276     I myOptimIter,
277     I mythid
278     & )
279    
280     C !DESCRIPTION: \bv
281     c ==================================================================
282     c SUBROUTINE active_write_rl
283     c ==================================================================
284     c o Write an active _RL variable to a file.
285     c started: Christian Eckert eckert@mit.edu Jan-1999
286     c ==================================================================
287     c SUBROUTINE active_write_rl
288     c ==================================================================
289     C \ev
290    
291     C !USES:
292     implicit none
293    
294     c == global variables ==
295     #include "EEPARAMS.h"
296     #include "SIZE.h"
297     #include "PARAMS.h"
298    
299     C !INPUT/OUTPUT PARAMETERS:
300     c == routine arguments ==
301     c active_var_file: filename
302     c active_var: array
303     c irec: record number
304     c myOptimIter: number of optimization iteration (default: 0)
305     c mythid: thread number for this instance
306     c doglobalread: flag for global or local read/write
307     c (default: .false.)
308     c lAdInit: initialisation of corresponding adjoint
309     c variable and write to active file
310     c mynr: vertical array dimension
311     c theSimulationMode: forward mode or reverse mode simulation
312     character*(*) active_var_file
313     integer mynr
314     logical globalfile
315     integer irec
316     integer theSimulationMode
317     integer myOptimIter
318     integer mythid
319     _RL active_var(1-olx:snx+olx,1-oly:sny+oly,mynr,nsx,nsy)
320    
321     C !LOCAL VARIABLES:
322     c == local variables ==
323     integer i,j,k
324     integer bi,bj
325     _RL active_data_t(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
326     integer oldprec
327     integer prec
328    
329     c == end of interface ==
330     CEOP
331    
332     c force 64-bit io
333     oldPrec = readBinaryPrec
334     readBinaryPrec = precFloat64
335     prec = precFloat64
336    
337     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
338     c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<<
339     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
340    
341     if (theSimulationMode .eq. FORWARD_SIMULATION) then
342    
343     _BEGIN_MASTER( mythid )
344    
345     call mdswritefield(
346     & active_var_file,
347     & prec,
348     & globalfile,
349     & 'RL',
350     & mynr,
351     & active_var,
352     & irec,
353     & myOptimIter,
354     & mythid )
355    
356     _END_MASTER( mythid )
357    
358     endif
359    
360     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
361     c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<<
362     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
363    
364     if (theSimulationMode .eq. REVERSE_SIMULATION) then
365    
366     _BEGIN_MASTER( mythid )
367    
368     do k=1,mynr
369     c Read data from file layer by layer.
370     call mdsreadfield(
371     & active_var_file,
372     & prec,
373     & 'RL',
374     & 1,
375     & active_data_t,
376     & (irec-1)*mynr+k,
377     & mythid )
378    
379     c Add active_var from appropriate location to data.
380     do bj = 1,nsy
381     do bi = 1,nsx
382     do j=1,sny
383     do i=1,snx
384     active_var(i,j,k,bi,bj) =
385     & active_var(i,j,k,bi,bj) +
386     & active_data_t(i,j,bi,bj)
387     active_data_t(i,j,bi,bj) = 0. _d 0
388     enddo
389     enddo
390     enddo
391     enddo
392     call mdswritefield(
393     & active_var_file,
394     & prec,
395     & globalfile,
396     & 'RL',
397     & 1,
398     & active_data_t,
399     & (irec-1)*mynr+k,
400     & myOptimIter,
401     & mythid )
402     enddo
403    
404     _END_MASTER( mythid )
405    
406     endif
407    
408 heimbach 1.4 c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
409     c >>>>>>>>>>>>>>>>>>> TANGENT RUN <<<<<<<<<<<<<<<<<<<
410     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
411    
412     if (theSimulationMode .eq. TANGENT_SIMULATION) then
413    
414     _BEGIN_MASTER( mythid )
415    
416     call mdswritefield(
417     & active_var_file,
418     & prec,
419     & globalfile,
420     & 'RL',
421     & mynr,
422     & active_var,
423     & irec,
424     & myOptimIter,
425     & mythid )
426    
427     _END_MASTER( mythid )
428    
429     endif
430    
431 heimbach 1.3 c Reset default io precision.
432     readBinaryPrec = oldPrec
433    
434     _BARRIER
435    
436     return
437     end
438    
439    
440     subroutine active_read_tile_rl(
441     I active_var_file,
442     O active_var,
443     I globalfile,
444     I lAdInit,
445     I irec,
446     I mynr,
447     I itile,
448     I jtile,
449     I theSimulationMode,
450     I myOptimIter,
451     I mythid
452     & )
453    
454     c ==================================================================
455     c SUBROUTINE active_read_tile_rl
456     c ==================================================================
457     c
458     c o Read an active variable from file.
459     c
460     c The variable *globalfile* can be used as a switch, which allows
461     c to read from a global file. The adjoint files are, however, always
462     c treated as tiled files.
463     c
464     c started: Christian Eckert eckert@mit.edu Jan-1999
465     c
466     c changed: Christian Eckert eckert@mit.edu 23-Nov-1999
467     c
468     c - interchanged the i,j loops
469     c
470     c Christian Eckert eckert@mit.edu 11-Feb-2000
471     c
472     c - Restructured the code in order to create a package
473     c for the MITgcmUV.
474     c
475     c Christian Eckert eckert@mit.edu 31-Mar-2000
476     c
477     c - Added BEGIN/END_MASTER and BARRIER for multiple
478     c threads and synchronisation, respectively. Added
479     c logical lAdInit.
480     c
481     c changed: Christian Eckert eckert@mit.edu 24-Apr-2000
482     c
483     c - Added routines that do active writes on tiles
484     c instead of a whole thread.
485     c
486     c ==================================================================
487     c SUBROUTINE active_read_tile_rl
488     c ==================================================================
489    
490     implicit none
491    
492     c == global variables ==
493    
494     #include "EEPARAMS.h"
495     #include "SIZE.h"
496     #include "PARAMS.h"
497    
498     c == routine arguments ==
499    
500     character*(*) active_var_file
501    
502     logical globalfile
503     logical lAdInit
504     integer irec
505     integer mynr
506     integer itile
507     integer jtile
508     integer theSimulationMode
509     integer myOptimIter
510     integer mythid
511     _RL active_var((snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy)
512    
513     c == local variables ==
514    
515     character*(2) adpref
516     character*(80) adfname
517    
518     integer bi,bj
519     integer i,j,k
520     integer index_t
521     integer index_var
522     integer oldprec
523     integer prec
524     integer il
525     integer ilnblnk
526    
527     logical writeglobalfile
528    
529     _RL active_data_t((snx+2*olx)*(sny+2*oly))
530    
531     c == functions ==
532    
533     external ilnblnk
534    
535     c == end of interface ==
536    
537     c force 64-bit io
538     oldPrec = readBinaryPrec
539     readBinaryPrec = precFloat64
540     prec = precFloat64
541    
542     write(adfname(1:80),'(80a)') ' '
543     adpref = 'ad'
544     il = ilnblnk( active_var_file )
545    
546     write(adfname(1:2) ,'(a)') adpref
547     write(adfname(3:il+2),'(a)') active_var_file(1:il)
548    
549     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
550     c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<<
551     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
552    
553     if (theSimulationMode .eq. FORWARD_SIMULATION) then
554    
555     _BEGIN_MASTER( mythid )
556    
557     c Read the active variable from file.
558     if ( .not. globalfile ) then
559     bj = jtile
560     bi = itile
561     c Read the layers individually.
562     do k=1,mynr
563     call mdsreadvector( active_var_file, prec, 'RL',
564     & (snx+2*olx)*(sny+2*oly),
565     & active_data_t, bi, bj,
566     & (irec-1)*mynr+k, mythid )
567    
568     c Copy data to appropriate location in active_var.
569     index_t = 0
570     do j=1-oly,sny+oly
571     do i=1-olx,snx+olx
572     index_t = index_t + 1
573     index_var = 1 +
574     & (i-1+olx) +
575     & (j-1+oly)*(snx+2*olx) +
576     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
577     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
578     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
579    
580     active_var(index_var) = active_data_t(index_t)
581     enddo
582     enddo
583     enddo
584    
585     else if ( globalfile ) then
586     do i=1,(snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy
587     active_var(i) = 0. _d 0
588     enddo
589     call mdsreadfield( active_var_file, prec, 'RL',
590     & mynr, active_var, irec, mythid )
591     endif
592    
593     if (lAdInit) then
594     c Initialise the corresponding adjoint variable on the
595     c adjoint variable's file.
596    
597     writeglobalfile = .false.
598     do i=1,(snx+2*olx)*(sny+2*oly)
599     active_data_t(i) = 0. _d 0
600     enddo
601     bj = jtile
602     bi = itile
603     do k = 1,mynr
604     call mdswritevector( adfname, prec, writeglobalfile,
605     & 'RL', (snx+2*olx)*(sny+2*oly),
606     & active_data_t, bi, bj,
607     & (irec-1)*mynr+k, myOptimIter,
608     & mythid )
609     enddo
610     endif
611    
612     _END_MASTER( mythid )
613    
614     endif
615    
616     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
617     c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<<
618     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
619    
620     if (theSimulationMode .eq. REVERSE_SIMULATION) then
621    
622     _BEGIN_MASTER( mythid )
623    
624     writeglobalfile = .false.
625     bj = jtile
626     bi = itile
627     do k=1,mynr
628     c Read data from file layer by layer.
629     call mdsreadvector( active_var_file, prec, 'RL',
630     & (snx+2*olx)*(sny+2*oly),
631     & active_data_t, bi, bj,
632     & (irec-1)*mynr+k, mythid )
633    
634     c Add active_var from appropriate location to data.
635     index_t = 0
636     do j=1-oly,sny+oly
637     do i=1-olx,snx+olx
638     index_t = index_t + 1
639     index_var = 1 +
640     & (i-1+olx) +
641     & (j-1+oly)*(snx+2*olx) +
642     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
643     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
644     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
645    
646     active_data_t(index_t) = active_data_t(index_t) +
647     & active_var(index_var)
648     enddo
649     enddo
650    
651     c Store the result on disk.
652     call mdswritevector( active_var_file, prec,
653     & writeglobalfile, 'RL',
654     & (snx+2*olx)*(sny+2*oly),
655     & active_data_t, bi, bj,
656     & (irec-1)*mynr+k, myOptimIter,
657     & mythid )
658     enddo
659    
660    
661     c Set active_var to zero. Use the standard FORTRAN index
662     c mapping.
663     bj = jtile
664     bi = itile
665     do k = 1,nr
666     do j = 1-oly,sny+oly
667     do i = 1-olx,snx+olx
668     index_var = 1 +
669     & (i-1+olx) +
670     & (j-1+oly)*(snx+2*olx) +
671     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
672     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
673     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
674    
675     active_var(index_var) = 0. _d 0
676     enddo
677     enddo
678     enddo
679    
680     _END_MASTER( mythid )
681    
682     endif
683    
684 heimbach 1.4 c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
685     c >>>>>>>>>>>>>>>>>>> TANGENT RUN <<<<<<<<<<<<<<<<<<<
686     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
687    
688     if (theSimulationMode .eq. TANGENT_SIMULATION) then
689    
690     _BEGIN_MASTER( mythid )
691    
692     c Read the active variable from file.
693     if ( .not. globalfile ) then
694     bj = jtile
695     bi = itile
696     c Read the layers individually.
697     do k=1,mynr
698     call mdsreadvector( active_var_file, prec, 'RL',
699     & (snx+2*olx)*(sny+2*oly),
700     & active_data_t, bi, bj,
701     & (irec-1)*mynr+k, mythid )
702    
703     c Copy data to appropriate location in active_var.
704     index_t = 0
705     do j=1-oly,sny+oly
706     do i=1-olx,snx+olx
707     index_t = index_t + 1
708     index_var = 1 +
709     & (i-1+olx) +
710     & (j-1+oly)*(snx+2*olx) +
711     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
712     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
713     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
714    
715     active_var(index_var) = active_data_t(index_t)
716     enddo
717     enddo
718     enddo
719    
720     else if ( globalfile ) then
721     do i=1,(snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy
722     active_var(i) = 0. _d 0
723     enddo
724     call mdsreadfield( active_var_file, prec, 'RL',
725     & mynr, active_var, irec, mythid )
726     endif
727    
728     _END_MASTER( mythid )
729    
730     endif
731    
732 heimbach 1.3 c Reset default io precision.
733     readBinaryPrec = oldPrec
734    
735     _BARRIER
736    
737     return
738     end
739    
740    
741     subroutine active_write_tile_rl(
742     I active_var_file,
743     I active_var,
744     I globalfile,
745     I irec,
746     I mynr,
747     I itile,
748     I jtile,
749     I theSimulationMode,
750     I myOptimIter,
751     I mythid
752     & )
753    
754     c ==================================================================
755     c SUBROUTINE active_write_tile_rl
756     c ==================================================================
757     c
758     c o Write an active variable to a file.
759     c
760     c started: Christian Eckert eckert@mit.edu Jan-1999
761     c
762     c changed: Christian Eckert eckert@mit.edu 23-Nov-1999
763     c
764     c - interchanged the i,j loops
765     c
766     c Christian Eckert eckert@mit.edu 11-Feb-2000
767     c
768     c - Restructured the code in order to create a package
769     c for the MITgcmUV.
770     c
771     c Christian Eckert eckert@mit.edu 31-Mar-2000
772     c
773     c - Added BEGIN/END_MASTER and BARRIER for multiple
774     c threads and synchronisation, respectively.c
775     c
776     c changed: Christian Eckert eckert@mit.edu 24-Apr-2000
777     c
778     c - Added routines that do active writes on tiles
779     c instead of a whole thread.
780     c
781     c ==================================================================
782     c SUBROUTINE active_write_tile_rl
783     c ==================================================================
784    
785     implicit none
786    
787     c == global variables ==
788    
789     #include "EEPARAMS.h"
790     #include "SIZE.h"
791     #include "PARAMS.h"
792    
793     c == routine arguments ==
794    
795     character*(*) active_var_file
796    
797     integer mynr
798     logical globalfile
799     integer irec
800     integer itile
801     integer jtile
802     integer theSimulationMode
803     integer myOptimIter
804     integer mythid
805     _RL active_var((snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy)
806    
807     c == local variables ==
808    
809     integer i,j,k
810     integer bi,bj
811     _RL active_data_t((snx+2*olx)*(sny+2*oly))
812     integer oldprec
813     integer prec
814     integer index_t
815     integer index_var
816    
817     c == end of interface ==
818    
819     c force 64-bit io
820     oldPrec = readBinaryPrec
821     readBinaryPrec = precFloat64
822     prec = precFloat64
823    
824     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
825     c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<<
826     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
827    
828     if (theSimulationMode .eq. FORWARD_SIMULATION) then
829    
830     _BEGIN_MASTER( mythid )
831    
832     bj = jtile
833     bi = itile
834     do k=1,mynr
835     index_t = 0
836     do j=1-oly,sny+oly
837     do i=1-olx,snx+olx
838     index_t = index_t + 1
839     index_var = 1 +
840     & (i-1+olx) +
841     & (j-1+oly)*(snx+2*olx) +
842     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
843     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
844     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
845    
846     active_data_t(index_t) = active_var(index_var)
847     enddo
848     enddo
849     call mdswritevector( active_var_file, prec, globalfile,
850     & 'RL', (snx+2*olx)*(sny+2*oly),
851     & active_data_t, bi, bj,
852     & (irec-1)*mynr+k, myOptimIter,
853     & mythid )
854    
855     enddo
856    
857     _END_MASTER( mythid )
858    
859     endif
860    
861     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
862     c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<<
863     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
864    
865     if (theSimulationMode .eq. REVERSE_SIMULATION) then
866    
867     _BEGIN_MASTER( mythid )
868    
869     bj = jtile
870     bi = itile
871     do k=1,mynr
872     c Read data from file layer by layer.
873     call mdsreadvector( active_var_file, prec, 'RL',
874     & (snx+2*olx)*(sny+2*oly),
875     & active_data_t, bi, bj,
876     & (irec-1)*mynr+k, mythid )
877    
878     c Add active_var from appropriate location to data.
879     index_t = 0
880     do j=1-oly,sny+oly
881     do i=1-olx,snx+olx
882     index_t = index_t + 1
883     index_var = 1 +
884     & (i-1+olx) +
885     & (j-1+oly)*(snx+2*olx) +
886     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
887     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
888     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
889    
890     active_var(index_var) = active_var(index_var) +
891     & active_data_t(index_t)
892     active_data_t(index_t) = 0. _d 0
893     enddo
894     enddo
895     call mdswritevector( active_var_file, prec, globalfile,
896     & 'RL', (snx+2*olx)*(sny+2*oly),
897     & active_data_t, bi, bj,
898     & (irec-1)*mynr+k, myOptimIter,
899     & mythid )
900 heimbach 1.4 enddo
901    
902     _END_MASTER( mythid )
903    
904     endif
905    
906     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
907     c >>>>>>>>>>>>>>>>>>> TANGENT RUN <<<<<<<<<<<<<<<<<<<
908     c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<<
909    
910     if (theSimulationMode .eq. TANGENT_SIMULATION) then
911    
912     _BEGIN_MASTER( mythid )
913    
914     bj = jtile
915     bi = itile
916     do k=1,mynr
917     index_t = 0
918     do j=1-oly,sny+oly
919     do i=1-olx,snx+olx
920     index_t = index_t + 1
921     index_var = 1 +
922     & (i-1+olx) +
923     & (j-1+oly)*(snx+2*olx) +
924     & ( k-1 )*(snx+2*olx)*(sny+2*oly) +
925     & ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr +
926     & ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx
927    
928     active_data_t(index_t) = active_var(index_var)
929     enddo
930     enddo
931     call mdswritevector( active_var_file, prec, globalfile,
932     & 'RL', (snx+2*olx)*(sny+2*oly),
933     & active_data_t, bi, bj,
934     & (irec-1)*mynr+k, myOptimIter,
935     & mythid )
936    
937 heimbach 1.3 enddo
938    
939     _END_MASTER( mythid )
940    
941     endif
942    
943     c Reset default io precision.
944     readBinaryPrec = oldPrec
945    
946     _BARRIER
947    
948     return
949     end

  ViewVC Help
Powered by ViewVC 1.1.22