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 - 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 |
c Reset default io precision. |
236 |
readBinaryPrec = oldPrec |
237 |
|
238 |
_BARRIER |
239 |
|
240 |
return |
241 |
end |
242 |
|
243 |
|
244 |
CBOP |
245 |
C !ROUTINE: active_write_rl |
246 |
C !INTERFACE: |
247 |
subroutine active_write_rl( |
248 |
I active_var_file, |
249 |
I active_var, |
250 |
I globalfile, |
251 |
I irec, |
252 |
I mynr, |
253 |
I theSimulationMode, |
254 |
I myOptimIter, |
255 |
I mythid |
256 |
& ) |
257 |
|
258 |
C !DESCRIPTION: \bv |
259 |
c ================================================================== |
260 |
c SUBROUTINE active_write_rl |
261 |
c ================================================================== |
262 |
c o Write an active _RL variable to a file. |
263 |
c started: Christian Eckert eckert@mit.edu Jan-1999 |
264 |
c ================================================================== |
265 |
c SUBROUTINE active_write_rl |
266 |
c ================================================================== |
267 |
C \ev |
268 |
|
269 |
C !USES: |
270 |
implicit none |
271 |
|
272 |
c == global variables == |
273 |
#include "EEPARAMS.h" |
274 |
#include "SIZE.h" |
275 |
#include "PARAMS.h" |
276 |
|
277 |
C !INPUT/OUTPUT PARAMETERS: |
278 |
c == routine arguments == |
279 |
c active_var_file: filename |
280 |
c active_var: array |
281 |
c irec: record number |
282 |
c myOptimIter: number of optimization iteration (default: 0) |
283 |
c mythid: thread number for this instance |
284 |
c doglobalread: flag for global or local read/write |
285 |
c (default: .false.) |
286 |
c lAdInit: initialisation of corresponding adjoint |
287 |
c variable and write to active file |
288 |
c mynr: vertical array dimension |
289 |
c theSimulationMode: forward mode or reverse mode simulation |
290 |
character*(*) active_var_file |
291 |
integer mynr |
292 |
logical globalfile |
293 |
integer irec |
294 |
integer theSimulationMode |
295 |
integer myOptimIter |
296 |
integer mythid |
297 |
_RL active_var(1-olx:snx+olx,1-oly:sny+oly,mynr,nsx,nsy) |
298 |
|
299 |
C !LOCAL VARIABLES: |
300 |
c == local variables == |
301 |
integer i,j,k |
302 |
integer bi,bj |
303 |
_RL active_data_t(1-olx:snx+olx,1-oly:sny+oly,nsx,nsy) |
304 |
integer oldprec |
305 |
integer prec |
306 |
|
307 |
c == end of interface == |
308 |
CEOP |
309 |
|
310 |
c force 64-bit io |
311 |
oldPrec = readBinaryPrec |
312 |
readBinaryPrec = precFloat64 |
313 |
prec = precFloat64 |
314 |
|
315 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
316 |
c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<< |
317 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
318 |
|
319 |
if (theSimulationMode .eq. FORWARD_SIMULATION) then |
320 |
|
321 |
_BEGIN_MASTER( mythid ) |
322 |
|
323 |
call mdswritefield( |
324 |
& active_var_file, |
325 |
& prec, |
326 |
& globalfile, |
327 |
& 'RL', |
328 |
& mynr, |
329 |
& active_var, |
330 |
& irec, |
331 |
& myOptimIter, |
332 |
& mythid ) |
333 |
|
334 |
_END_MASTER( mythid ) |
335 |
|
336 |
endif |
337 |
|
338 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
339 |
c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<< |
340 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
341 |
|
342 |
if (theSimulationMode .eq. REVERSE_SIMULATION) then |
343 |
|
344 |
_BEGIN_MASTER( mythid ) |
345 |
|
346 |
do k=1,mynr |
347 |
c Read data from file layer by layer. |
348 |
call mdsreadfield( |
349 |
& active_var_file, |
350 |
& prec, |
351 |
& 'RL', |
352 |
& 1, |
353 |
& active_data_t, |
354 |
& (irec-1)*mynr+k, |
355 |
& mythid ) |
356 |
|
357 |
c Add active_var from appropriate location to data. |
358 |
do bj = 1,nsy |
359 |
do bi = 1,nsx |
360 |
do j=1,sny |
361 |
do i=1,snx |
362 |
active_var(i,j,k,bi,bj) = |
363 |
& active_var(i,j,k,bi,bj) + |
364 |
& active_data_t(i,j,bi,bj) |
365 |
active_data_t(i,j,bi,bj) = 0. _d 0 |
366 |
enddo |
367 |
enddo |
368 |
enddo |
369 |
enddo |
370 |
call mdswritefield( |
371 |
& active_var_file, |
372 |
& prec, |
373 |
& globalfile, |
374 |
& 'RL', |
375 |
& 1, |
376 |
& active_data_t, |
377 |
& (irec-1)*mynr+k, |
378 |
& myOptimIter, |
379 |
& mythid ) |
380 |
enddo |
381 |
|
382 |
_END_MASTER( mythid ) |
383 |
|
384 |
endif |
385 |
|
386 |
c Reset default io precision. |
387 |
readBinaryPrec = oldPrec |
388 |
|
389 |
_BARRIER |
390 |
|
391 |
return |
392 |
end |
393 |
|
394 |
|
395 |
subroutine active_read_tile_rl( |
396 |
I active_var_file, |
397 |
O active_var, |
398 |
I globalfile, |
399 |
I lAdInit, |
400 |
I irec, |
401 |
I mynr, |
402 |
I itile, |
403 |
I jtile, |
404 |
I theSimulationMode, |
405 |
I myOptimIter, |
406 |
I mythid |
407 |
& ) |
408 |
|
409 |
c ================================================================== |
410 |
c SUBROUTINE active_read_tile_rl |
411 |
c ================================================================== |
412 |
c |
413 |
c o Read an active variable from file. |
414 |
c |
415 |
c The variable *globalfile* can be used as a switch, which allows |
416 |
c to read from a global file. The adjoint files are, however, always |
417 |
c treated as tiled files. |
418 |
c |
419 |
c started: Christian Eckert eckert@mit.edu Jan-1999 |
420 |
c |
421 |
c changed: Christian Eckert eckert@mit.edu 23-Nov-1999 |
422 |
c |
423 |
c - interchanged the i,j loops |
424 |
c |
425 |
c Christian Eckert eckert@mit.edu 11-Feb-2000 |
426 |
c |
427 |
c - Restructured the code in order to create a package |
428 |
c for the MITgcmUV. |
429 |
c |
430 |
c Christian Eckert eckert@mit.edu 31-Mar-2000 |
431 |
c |
432 |
c - Added BEGIN/END_MASTER and BARRIER for multiple |
433 |
c threads and synchronisation, respectively. Added |
434 |
c logical lAdInit. |
435 |
c |
436 |
c changed: Christian Eckert eckert@mit.edu 24-Apr-2000 |
437 |
c |
438 |
c - Added routines that do active writes on tiles |
439 |
c instead of a whole thread. |
440 |
c |
441 |
c ================================================================== |
442 |
c SUBROUTINE active_read_tile_rl |
443 |
c ================================================================== |
444 |
|
445 |
implicit none |
446 |
|
447 |
c == global variables == |
448 |
|
449 |
#include "EEPARAMS.h" |
450 |
#include "SIZE.h" |
451 |
#include "PARAMS.h" |
452 |
|
453 |
c == routine arguments == |
454 |
|
455 |
character*(*) active_var_file |
456 |
|
457 |
logical globalfile |
458 |
logical lAdInit |
459 |
integer irec |
460 |
integer mynr |
461 |
integer itile |
462 |
integer jtile |
463 |
integer theSimulationMode |
464 |
integer myOptimIter |
465 |
integer mythid |
466 |
_RL active_var((snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy) |
467 |
|
468 |
c == local variables == |
469 |
|
470 |
character*(2) adpref |
471 |
character*(80) adfname |
472 |
|
473 |
integer bi,bj |
474 |
integer i,j,k |
475 |
integer index_t |
476 |
integer index_var |
477 |
integer oldprec |
478 |
integer prec |
479 |
integer il |
480 |
integer ilnblnk |
481 |
|
482 |
logical writeglobalfile |
483 |
|
484 |
_RL active_data_t((snx+2*olx)*(sny+2*oly)) |
485 |
|
486 |
c == functions == |
487 |
|
488 |
external ilnblnk |
489 |
|
490 |
c == end of interface == |
491 |
|
492 |
c force 64-bit io |
493 |
oldPrec = readBinaryPrec |
494 |
readBinaryPrec = precFloat64 |
495 |
prec = precFloat64 |
496 |
|
497 |
write(adfname(1:80),'(80a)') ' ' |
498 |
adpref = 'ad' |
499 |
il = ilnblnk( active_var_file ) |
500 |
|
501 |
write(adfname(1:2) ,'(a)') adpref |
502 |
write(adfname(3:il+2),'(a)') active_var_file(1:il) |
503 |
|
504 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
505 |
c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<< |
506 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
507 |
|
508 |
if (theSimulationMode .eq. FORWARD_SIMULATION) then |
509 |
|
510 |
_BEGIN_MASTER( mythid ) |
511 |
|
512 |
c Read the active variable from file. |
513 |
if ( .not. globalfile ) then |
514 |
bj = jtile |
515 |
bi = itile |
516 |
c Read the layers individually. |
517 |
do k=1,mynr |
518 |
call mdsreadvector( active_var_file, prec, 'RL', |
519 |
& (snx+2*olx)*(sny+2*oly), |
520 |
& active_data_t, bi, bj, |
521 |
& (irec-1)*mynr+k, mythid ) |
522 |
|
523 |
c Copy data to appropriate location in active_var. |
524 |
index_t = 0 |
525 |
do j=1-oly,sny+oly |
526 |
do i=1-olx,snx+olx |
527 |
index_t = index_t + 1 |
528 |
index_var = 1 + |
529 |
& (i-1+olx) + |
530 |
& (j-1+oly)*(snx+2*olx) + |
531 |
& ( k-1 )*(snx+2*olx)*(sny+2*oly) + |
532 |
& ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr + |
533 |
& ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx |
534 |
|
535 |
active_var(index_var) = active_data_t(index_t) |
536 |
enddo |
537 |
enddo |
538 |
enddo |
539 |
|
540 |
else if ( globalfile ) then |
541 |
do i=1,(snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy |
542 |
active_var(i) = 0. _d 0 |
543 |
enddo |
544 |
call mdsreadfield( active_var_file, prec, 'RL', |
545 |
& mynr, active_var, irec, mythid ) |
546 |
endif |
547 |
|
548 |
if (lAdInit) then |
549 |
c Initialise the corresponding adjoint variable on the |
550 |
c adjoint variable's file. |
551 |
|
552 |
writeglobalfile = .false. |
553 |
do i=1,(snx+2*olx)*(sny+2*oly) |
554 |
active_data_t(i) = 0. _d 0 |
555 |
enddo |
556 |
bj = jtile |
557 |
bi = itile |
558 |
do k = 1,mynr |
559 |
call mdswritevector( adfname, prec, writeglobalfile, |
560 |
& 'RL', (snx+2*olx)*(sny+2*oly), |
561 |
& active_data_t, bi, bj, |
562 |
& (irec-1)*mynr+k, myOptimIter, |
563 |
& mythid ) |
564 |
enddo |
565 |
endif |
566 |
|
567 |
_END_MASTER( mythid ) |
568 |
|
569 |
endif |
570 |
|
571 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
572 |
c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<< |
573 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
574 |
|
575 |
if (theSimulationMode .eq. REVERSE_SIMULATION) then |
576 |
|
577 |
_BEGIN_MASTER( mythid ) |
578 |
|
579 |
writeglobalfile = .false. |
580 |
bj = jtile |
581 |
bi = itile |
582 |
do k=1,mynr |
583 |
c Read data from file layer by layer. |
584 |
call mdsreadvector( active_var_file, prec, 'RL', |
585 |
& (snx+2*olx)*(sny+2*oly), |
586 |
& active_data_t, bi, bj, |
587 |
& (irec-1)*mynr+k, mythid ) |
588 |
|
589 |
c Add active_var from appropriate location to data. |
590 |
index_t = 0 |
591 |
do j=1-oly,sny+oly |
592 |
do i=1-olx,snx+olx |
593 |
index_t = index_t + 1 |
594 |
index_var = 1 + |
595 |
& (i-1+olx) + |
596 |
& (j-1+oly)*(snx+2*olx) + |
597 |
& ( k-1 )*(snx+2*olx)*(sny+2*oly) + |
598 |
& ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr + |
599 |
& ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx |
600 |
|
601 |
active_data_t(index_t) = active_data_t(index_t) + |
602 |
& active_var(index_var) |
603 |
enddo |
604 |
enddo |
605 |
|
606 |
c Store the result on disk. |
607 |
call mdswritevector( active_var_file, prec, |
608 |
& writeglobalfile, 'RL', |
609 |
& (snx+2*olx)*(sny+2*oly), |
610 |
& active_data_t, bi, bj, |
611 |
& (irec-1)*mynr+k, myOptimIter, |
612 |
& mythid ) |
613 |
enddo |
614 |
|
615 |
|
616 |
c Set active_var to zero. Use the standard FORTRAN index |
617 |
c mapping. |
618 |
bj = jtile |
619 |
bi = itile |
620 |
do k = 1,nr |
621 |
do j = 1-oly,sny+oly |
622 |
do i = 1-olx,snx+olx |
623 |
index_var = 1 + |
624 |
& (i-1+olx) + |
625 |
& (j-1+oly)*(snx+2*olx) + |
626 |
& ( k-1 )*(snx+2*olx)*(sny+2*oly) + |
627 |
& ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr + |
628 |
& ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx |
629 |
|
630 |
active_var(index_var) = 0. _d 0 |
631 |
enddo |
632 |
enddo |
633 |
enddo |
634 |
|
635 |
_END_MASTER( mythid ) |
636 |
|
637 |
endif |
638 |
|
639 |
c Reset default io precision. |
640 |
readBinaryPrec = oldPrec |
641 |
|
642 |
_BARRIER |
643 |
|
644 |
return |
645 |
end |
646 |
|
647 |
|
648 |
subroutine active_write_tile_rl( |
649 |
I active_var_file, |
650 |
I active_var, |
651 |
I globalfile, |
652 |
I irec, |
653 |
I mynr, |
654 |
I itile, |
655 |
I jtile, |
656 |
I theSimulationMode, |
657 |
I myOptimIter, |
658 |
I mythid |
659 |
& ) |
660 |
|
661 |
c ================================================================== |
662 |
c SUBROUTINE active_write_tile_rl |
663 |
c ================================================================== |
664 |
c |
665 |
c o Write an active variable to a file. |
666 |
c |
667 |
c started: Christian Eckert eckert@mit.edu Jan-1999 |
668 |
c |
669 |
c changed: Christian Eckert eckert@mit.edu 23-Nov-1999 |
670 |
c |
671 |
c - interchanged the i,j loops |
672 |
c |
673 |
c Christian Eckert eckert@mit.edu 11-Feb-2000 |
674 |
c |
675 |
c - Restructured the code in order to create a package |
676 |
c for the MITgcmUV. |
677 |
c |
678 |
c Christian Eckert eckert@mit.edu 31-Mar-2000 |
679 |
c |
680 |
c - Added BEGIN/END_MASTER and BARRIER for multiple |
681 |
c threads and synchronisation, respectively.c |
682 |
c |
683 |
c changed: Christian Eckert eckert@mit.edu 24-Apr-2000 |
684 |
c |
685 |
c - Added routines that do active writes on tiles |
686 |
c instead of a whole thread. |
687 |
c |
688 |
c ================================================================== |
689 |
c SUBROUTINE active_write_tile_rl |
690 |
c ================================================================== |
691 |
|
692 |
implicit none |
693 |
|
694 |
c == global variables == |
695 |
|
696 |
#include "EEPARAMS.h" |
697 |
#include "SIZE.h" |
698 |
#include "PARAMS.h" |
699 |
|
700 |
c == routine arguments == |
701 |
|
702 |
character*(*) active_var_file |
703 |
|
704 |
integer mynr |
705 |
logical globalfile |
706 |
integer irec |
707 |
integer itile |
708 |
integer jtile |
709 |
integer theSimulationMode |
710 |
integer myOptimIter |
711 |
integer mythid |
712 |
_RL active_var((snx+2*olx)*(sny+2*oly)*mynr*nsx*nsy) |
713 |
|
714 |
c == local variables == |
715 |
|
716 |
integer i,j,k |
717 |
integer bi,bj |
718 |
_RL active_data_t((snx+2*olx)*(sny+2*oly)) |
719 |
integer oldprec |
720 |
integer prec |
721 |
integer index_t |
722 |
integer index_var |
723 |
|
724 |
c == end of interface == |
725 |
|
726 |
c force 64-bit io |
727 |
oldPrec = readBinaryPrec |
728 |
readBinaryPrec = precFloat64 |
729 |
prec = precFloat64 |
730 |
|
731 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
732 |
c >>>>>>>>>>>>>>>>>>> FORWARD RUN <<<<<<<<<<<<<<<<<<< |
733 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
734 |
|
735 |
if (theSimulationMode .eq. FORWARD_SIMULATION) then |
736 |
|
737 |
_BEGIN_MASTER( mythid ) |
738 |
|
739 |
bj = jtile |
740 |
bi = itile |
741 |
do k=1,mynr |
742 |
index_t = 0 |
743 |
do j=1-oly,sny+oly |
744 |
do i=1-olx,snx+olx |
745 |
index_t = index_t + 1 |
746 |
index_var = 1 + |
747 |
& (i-1+olx) + |
748 |
& (j-1+oly)*(snx+2*olx) + |
749 |
& ( k-1 )*(snx+2*olx)*(sny+2*oly) + |
750 |
& ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr + |
751 |
& ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx |
752 |
|
753 |
active_data_t(index_t) = active_var(index_var) |
754 |
enddo |
755 |
enddo |
756 |
call mdswritevector( active_var_file, prec, globalfile, |
757 |
& 'RL', (snx+2*olx)*(sny+2*oly), |
758 |
& active_data_t, bi, bj, |
759 |
& (irec-1)*mynr+k, myOptimIter, |
760 |
& mythid ) |
761 |
|
762 |
enddo |
763 |
|
764 |
_END_MASTER( mythid ) |
765 |
|
766 |
endif |
767 |
|
768 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
769 |
c >>>>>>>>>>>>>>>>>>> ADJOINT RUN <<<<<<<<<<<<<<<<<<< |
770 |
c >>>>>>>>>>>>>>>>>>> <<<<<<<<<<<<<<<<<<< |
771 |
|
772 |
if (theSimulationMode .eq. REVERSE_SIMULATION) then |
773 |
|
774 |
_BEGIN_MASTER( mythid ) |
775 |
|
776 |
bj = jtile |
777 |
bi = itile |
778 |
do k=1,mynr |
779 |
c Read data from file layer by layer. |
780 |
call mdsreadvector( active_var_file, prec, 'RL', |
781 |
& (snx+2*olx)*(sny+2*oly), |
782 |
& active_data_t, bi, bj, |
783 |
& (irec-1)*mynr+k, mythid ) |
784 |
|
785 |
c Add active_var from appropriate location to data. |
786 |
index_t = 0 |
787 |
do j=1-oly,sny+oly |
788 |
do i=1-olx,snx+olx |
789 |
index_t = index_t + 1 |
790 |
index_var = 1 + |
791 |
& (i-1+olx) + |
792 |
& (j-1+oly)*(snx+2*olx) + |
793 |
& ( k-1 )*(snx+2*olx)*(sny+2*oly) + |
794 |
& ( bi-1 )*(snx+2*olx)*(sny+2*oly)*nr + |
795 |
& ( bj-1 )*(snx+2*olx)*(sny+2*oly)*nr*nsx |
796 |
|
797 |
active_var(index_var) = active_var(index_var) + |
798 |
& active_data_t(index_t) |
799 |
active_data_t(index_t) = 0. _d 0 |
800 |
enddo |
801 |
enddo |
802 |
call mdswritevector( active_var_file, prec, globalfile, |
803 |
& 'RL', (snx+2*olx)*(sny+2*oly), |
804 |
& active_data_t, bi, bj, |
805 |
& (irec-1)*mynr+k, myOptimIter, |
806 |
& mythid ) |
807 |
enddo |
808 |
|
809 |
_END_MASTER( mythid ) |
810 |
|
811 |
endif |
812 |
|
813 |
c Reset default io precision. |
814 |
readBinaryPrec = oldPrec |
815 |
|
816 |
_BARRIER |
817 |
|
818 |
return |
819 |
end |