1 |
C $Header: /u/gcmpack/MITgcm/pkg/mdsio/mdsio_writevector.F,v 1.4 2005/07/20 20:01:46 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "MDSIO_OPTIONS.h" |
5 |
|
6 |
SUBROUTINE MDSWRITEVECTOR( |
7 |
I fName, |
8 |
I filePrec, |
9 |
I globalfile, |
10 |
I arrType, |
11 |
I narr, |
12 |
I arr, |
13 |
I bi, |
14 |
I bj, |
15 |
I irecord, |
16 |
I myIter, |
17 |
I myThid ) |
18 |
C Arguments: |
19 |
C |
20 |
C fName string base name for file to written |
21 |
C filePrec integer number of bits per word in file (32 or 64) |
22 |
C globalFile logical selects between writing a global or tiled file |
23 |
C arrType char(2) declaration of "arr": either "RS" or "RL" |
24 |
C narr integer size of third dimension: normally either 1 or Nr |
25 |
C arr RS/RL array to write, arr(narr) |
26 |
ce bi integer x tile index |
27 |
ce bj integer y tile index |
28 |
C irecord integer record number to read |
29 |
C myIter integer time step number |
30 |
C myThid integer thread identifier |
31 |
C |
32 |
C Created: 03/26/99 eckert@mit.edu |
33 |
C Modified: 03/29/99 adcroft@mit.edu + eckert@mit.edu |
34 |
C Fixed to work work with _RS and _RL declarations |
35 |
C Modified: 07/27/99 eckert@mit.edu |
36 |
C Customized for state estimation (--> active_file_control.F) |
37 |
C Changed: 05/31/00 heimbach@mit.edu |
38 |
C open(dUnit, ..., status='old', ... -> status='unknown' |
39 |
|
40 |
implicit none |
41 |
C Global variables / common blocks |
42 |
#include "SIZE.h" |
43 |
#include "EEPARAMS.h" |
44 |
#include "PARAMS.h" |
45 |
#include "EESUPPORT.h" |
46 |
|
47 |
C Routine arguments |
48 |
character*(*) fName |
49 |
integer filePrec |
50 |
logical globalfile |
51 |
character*(2) arrType |
52 |
integer narr |
53 |
Real arr(narr) |
54 |
integer irecord |
55 |
integer myIter |
56 |
integer myThid |
57 |
ce |
58 |
integer bi,bj |
59 |
ce |
60 |
|
61 |
C Functions |
62 |
integer ILNBLNK |
63 |
integer MDS_RECLEN |
64 |
C Local variables |
65 |
character*(128) dataFName,metaFName,pfName |
66 |
integer iG,jG,irec,dUnit,IL,pIL |
67 |
logical fileIsOpen |
68 |
integer dimList(3,3),ndims |
69 |
integer length_of_rec |
70 |
character*(max_len_mbuf) msgbuf |
71 |
|
72 |
cph( |
73 |
cph Deal with useSingleCpuIO |
74 |
cph Not implemented here for EXCH2 |
75 |
logical lprint |
76 |
INTEGER K,L |
77 |
INTEGER nNz |
78 |
INTEGER vec_size |
79 |
#ifdef ALLOW_USE_MPI |
80 |
INTEGER iG_IO,jG_IO,npe |
81 |
Real*4 xy_buffer_r4(narr*nPx*nPy) |
82 |
Real*8 xy_buffer_r8(narr*nPx*nPy) |
83 |
Real*8 global(narr*nPx*nPy) |
84 |
_RL local(narr) |
85 |
#endif |
86 |
cph) |
87 |
|
88 |
C ------------------------------------------------------------------ |
89 |
|
90 |
vec_size = narr*nPx*nPy |
91 |
nNz = 1 |
92 |
|
93 |
C Only do I/O if I am the master thread |
94 |
_BEGIN_MASTER( myThid ) |
95 |
|
96 |
C Record number must be >= 1 |
97 |
if (irecord .LT. 1) then |
98 |
write(msgbuf,'(a,i9.8)') |
99 |
& ' MDSWRITEVECTOR: argument irecord = ',irecord |
100 |
call print_message( msgbuf, standardmessageunit, |
101 |
& SQUEEZE_RIGHT , mythid) |
102 |
write(msgbuf,'(a)') |
103 |
& ' MDSWRITEVECTOR: invalid value for irecord' |
104 |
call print_error( msgbuf, mythid ) |
105 |
stop 'ABNORMAL END: S/R MDSWRITEVECTOR' |
106 |
endif |
107 |
|
108 |
C Assume nothing |
109 |
fileIsOpen = .FALSE. |
110 |
IL = ILNBLNK( fName ) |
111 |
pIL = ILNBLNK( mdsioLocalDir ) |
112 |
|
113 |
C Assign special directory |
114 |
if ( mdsioLocalDir .NE. ' ' ) then |
115 |
write(pFname(1:128),'(2a)') |
116 |
& mdsioLocalDir(1:pIL), fName(1:IL) |
117 |
else |
118 |
pFname= fName |
119 |
endif |
120 |
pIL=ILNBLNK( pfName ) |
121 |
|
122 |
C Assign a free unit number as the I/O channel for this routine |
123 |
call MDSFINDUNIT( dUnit, mythid ) |
124 |
|
125 |
#ifdef ALLOW_USE_MPI |
126 |
_END_MASTER( myThid ) |
127 |
C If option globalFile is desired but does not work or if |
128 |
C globalFile is too slow, then try using single-CPU I/O. |
129 |
if (useSingleCpuIO) then |
130 |
|
131 |
C Master thread of process 0, only, opens a global file |
132 |
_BEGIN_MASTER( myThid ) |
133 |
IF( mpiMyId .EQ. 0 ) THEN |
134 |
write(dataFname(1:128),'(2a)') fName(1:IL),'.data' |
135 |
length_of_rec=MDS_RECLEN(filePrec,vec_size,mythid) |
136 |
if (irecord .EQ. 1) then |
137 |
open( dUnit, file=dataFName, status=_NEW_STATUS, |
138 |
& access='direct', recl=length_of_rec ) |
139 |
else |
140 |
open( dUnit, file=dataFName, status=_OLD_STATUS, |
141 |
& access='direct', recl=length_of_rec ) |
142 |
endif |
143 |
ENDIF |
144 |
_END_MASTER( myThid ) |
145 |
|
146 |
C Gather array and write it to file, one vertical level at a time |
147 |
DO k=1,1 |
148 |
DO L=1,narr |
149 |
local(L) = arr(L) |
150 |
ENDDO |
151 |
cph( |
152 |
cph if ( irecord .EQ. 1 .AND. fName(1:IL) .EQ. |
153 |
cph & 'tapelev2_7_the_main_loop_theta.it0000' ) then |
154 |
cph lprint = .TRUE. |
155 |
cph else |
156 |
lprint = .FALSE. |
157 |
cph endif |
158 |
cph) |
159 |
CALL GATHER_VECTOR( lprint, narr, global, local, myThid ) |
160 |
_BEGIN_MASTER( myThid ) |
161 |
IF( mpiMyId .EQ. 0 ) THEN |
162 |
irec=irecord |
163 |
if (filePrec .eq. precFloat32) then |
164 |
#if defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) |
165 |
c |
166 |
#else /* defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) */ |
167 |
DO L=1,narr*nPx*nPy |
168 |
xy_buffer_r4(L) = global(L) |
169 |
ENDDO |
170 |
#endif /* defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) */ |
171 |
#ifdef _BYTESWAPIO |
172 |
call MDS_BYTESWAPR4( vec_size, xy_buffer_r4 ) |
173 |
#endif |
174 |
write(dUnit,rec=irec) xy_buffer_r4 |
175 |
elseif (filePrec .eq. precFloat64) then |
176 |
#if defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) |
177 |
c |
178 |
#else /* defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) */ |
179 |
DO L=1,narr*nPx*nPy |
180 |
xy_buffer_r8(L) = global(L) |
181 |
ENDDO |
182 |
#endif /* defined(ALLOW_EXCH2) && !defined(MISSING_TILE_IO) */ |
183 |
#ifdef _BYTESWAPIO |
184 |
call MDS_BYTESWAPR8( vec_size, xy_buffer_r8 ) |
185 |
#endif |
186 |
write(dUnit,rec=irec) xy_buffer_r8 |
187 |
else |
188 |
write(msgbuf,'(a)') |
189 |
& ' MDSWRITEFIELD: illegal value for filePrec' |
190 |
call print_error( msgbuf, mythid ) |
191 |
stop 'ABNORMAL END: S/R MDSWRITEFIELD' |
192 |
endif |
193 |
ENDIF |
194 |
_END_MASTER( myThid ) |
195 |
ENDDO |
196 |
|
197 |
C Close data-file and create meta-file |
198 |
_BEGIN_MASTER( myThid ) |
199 |
IF( mpiMyId .EQ. 0 ) THEN |
200 |
close( dUnit ) |
201 |
write(metaFName(1:128),'(2a)') fName(1:IL),'.meta' |
202 |
dimList(1,1)=vec_size |
203 |
dimList(2,1)=1 |
204 |
dimList(3,1)=vec_size |
205 |
dimList(1,2)=vec_size |
206 |
dimList(2,2)=1 |
207 |
dimList(3,2)=vec_size |
208 |
dimList(1,3)=1 |
209 |
dimList(2,3)=1 |
210 |
dimList(3,3)=1 |
211 |
ndims=1 |
212 |
cph if (nNz .EQ. 1) ndims=2 |
213 |
call MDSWRITEMETA( metaFName, dataFName, |
214 |
& filePrec, ndims, dimList, irecord, myIter, mythid ) |
215 |
ENDIF |
216 |
_END_MASTER( myThid ) |
217 |
C To be safe, make other processes wait for I/O completion |
218 |
_BARRIER |
219 |
|
220 |
elseif ( .NOT. useSingleCpuIO ) then |
221 |
_BEGIN_MASTER( myThid ) |
222 |
#endif /* ALLOW_USE_MPI */ |
223 |
|
224 |
C If we are writing to a global file then we open it here |
225 |
if (globalFile) then |
226 |
write(dataFname(1:128),'(2a)') fName(1:IL),'.data' |
227 |
if (irecord .EQ. 1) then |
228 |
length_of_rec = MDS_RECLEN( filePrec, narr, mythid ) |
229 |
open( dUnit, file=dataFName, status=_NEW_STATUS, |
230 |
& access='direct', recl=length_of_rec ) |
231 |
fileIsOpen=.TRUE. |
232 |
else |
233 |
length_of_rec = MDS_RECLEN( filePrec, narr, mythid ) |
234 |
open( dUnit, file=dataFName, status=_OLD_STATUS, |
235 |
& access='direct', recl=length_of_rec ) |
236 |
fileIsOpen=.TRUE. |
237 |
endif |
238 |
endif |
239 |
|
240 |
C Loop over all tiles |
241 |
ce do bj=1,nSy |
242 |
ce do bi=1,nSx |
243 |
C If we are writing to a tiled MDS file then we open each one here |
244 |
if (.NOT. globalFile) then |
245 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
246 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
247 |
write(dataFname(1:128),'(2a,i3.3,a,i3.3,a)') |
248 |
& pfName(1:pIL),'.',iG,'.',jG,'.data' |
249 |
if (irecord .EQ. 1) then |
250 |
length_of_rec = MDS_RECLEN( filePrec, narr, mythid ) |
251 |
open( dUnit, file=dataFName, status=_NEW_STATUS, |
252 |
& access='direct', recl=length_of_rec ) |
253 |
fileIsOpen=.TRUE. |
254 |
else |
255 |
length_of_rec = MDS_RECLEN( filePrec, narr, mythid ) |
256 |
open( dUnit, file=dataFName, status=_OLD_STATUS, |
257 |
& access='direct', recl=length_of_rec ) |
258 |
fileIsOpen=.TRUE. |
259 |
endif |
260 |
endif |
261 |
if (fileIsOpen) then |
262 |
if (globalFile) then |
263 |
iG = myXGlobalLo-1+(bi-1)*sNx |
264 |
jG = myYGlobalLo-1+(bj-1)*sNy |
265 |
irec = 1 + int(iG/sNx) + (jG/sNy)*nSx*nPx + |
266 |
& (irecord-1)*nSx*nPx*nSy*nPy |
267 |
else |
268 |
iG = 0 |
269 |
jG = 0 |
270 |
irec = irecord |
271 |
endif |
272 |
if (filePrec .eq. precFloat32) then |
273 |
call MDS_WRITE_RS_VEC( dUnit, irec, narr, arr, myThid ) |
274 |
elseif (filePrec .eq. precFloat64) then |
275 |
call MDS_WRITE_RL_VEC( dUnit, irec, narr, arr, myThid ) |
276 |
else |
277 |
write(msgbuf,'(a)') |
278 |
& ' MDSWRITEVECTOR: illegal value for filePrec' |
279 |
call print_error( msgbuf, mythid ) |
280 |
stop 'ABNORMAL END: S/R MDSWRITEVECTOR' |
281 |
endif |
282 |
else |
283 |
write(msgbuf,'(a)') |
284 |
& ' MDSWRITEVECTOR: I should never get to this point' |
285 |
call print_error( msgbuf, mythid ) |
286 |
stop 'ABNORMAL END: S/R MDSWRITEVECTOR' |
287 |
endif |
288 |
C If we were writing to a tiled MDS file then we close it here |
289 |
if (fileIsOpen .AND. (.NOT. globalFile)) then |
290 |
close( dUnit ) |
291 |
fileIsOpen = .FALSE. |
292 |
endif |
293 |
C Create meta-file for each tile file |
294 |
if (.NOT. globalFile) then |
295 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
296 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
297 |
write(metaFname(1:128),'(2a,i3.3,a,i3.3,a)') |
298 |
& pfName(1:pIL),'.',iG,'.',jG,'.meta' |
299 |
dimList(1,1) = nPx*nSx*narr |
300 |
dimList(2,1) = ((myXGlobalLo-1)/sNx + (bi-1))*narr + 1 |
301 |
dimList(3,1) = ((myXGlobalLo-1)/sNx + bi )*narr |
302 |
dimList(1,2) = nPy*nSy |
303 |
dimList(2,2) = (myYGlobalLo-1)/sNy + bj |
304 |
dimList(3,2) = (myYGlobalLo-1)/sNy + bj |
305 |
dimList(1,3) = 1 |
306 |
dimList(2,3) = 1 |
307 |
dimList(3,3) = 1 |
308 |
ndims=1 |
309 |
call MDSWRITEMETA( metaFName, dataFName, |
310 |
& filePrec, ndims, dimList, irecord, myIter, mythid ) |
311 |
endif |
312 |
C End of bi,bj loops |
313 |
ce enddo |
314 |
ce enddo |
315 |
|
316 |
C If global file was opened then close it |
317 |
if (fileIsOpen .AND. globalFile) then |
318 |
close( dUnit ) |
319 |
fileIsOpen = .FALSE. |
320 |
endif |
321 |
|
322 |
C Create meta-file for global file |
323 |
if (globalFile) then |
324 |
write(metaFName(1:128),'(2a)') fName(1:IL),'.meta' |
325 |
dimList(1,1) = nPx*nSx*narr |
326 |
dimList(2,1) = 1 |
327 |
dimList(3,1) = nPx*nSx*narr |
328 |
dimList(1,2) = nPy*nSy |
329 |
dimList(2,2) = 1 |
330 |
dimList(3,2) = nPy*nSy |
331 |
dimList(1,3) = 1 |
332 |
dimList(2,3) = 1 |
333 |
dimList(3,3) = 1 |
334 |
ndims=1 |
335 |
call MDSWRITEMETA( metaFName, dataFName, |
336 |
& filePrec, ndims, dimList, irecord, myIter, mythid ) |
337 |
endif |
338 |
|
339 |
_END_MASTER( myThid ) |
340 |
|
341 |
#ifdef ALLOW_USE_MPI |
342 |
C endif useSingleCpuIO |
343 |
endif |
344 |
#endif /* ALLOW_USE_MPI */ |
345 |
|
346 |
C ------------------------------------------------------------------ |
347 |
return |
348 |
end |