/[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.1 - (hide annotations) (download)
Sun Mar 25 22:33:53 2001 UTC (23 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint38, checkpoint40pre2, checkpoint40pre4, c37_adj, checkpoint39, checkpoint40pre5, checkpoint42, checkpoint40, checkpoint41
Modifications and additions to enable automatic differentiation.
Detailed info's in doc/notes_c37_adj.txt

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

  ViewVC Help
Powered by ViewVC 1.1.22