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

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

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


Revision 1.9 - (show annotations) (download)
Mon Oct 8 23:50:52 2007 UTC (16 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62a, checkpoint60, checkpoint61, checkpoint62, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59k, checkpoint62b, checkpoint61f, checkpoint61n, checkpoint59j, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.8: +59 -57 lines
add missing cvs $Header:$ or $Name:$

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

  ViewVC Help
Powered by ViewVC 1.1.22