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