1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "MDSIO_OPTIONS.h" |
5 |
|
6 |
SUBROUTINE MDS_WRITETILE( |
7 |
I fName, |
8 |
I filePrec, |
9 |
I globalFile, |
10 |
I arrType, |
11 |
I nNz, |
12 |
I arr, |
13 |
I bi, bj, |
14 |
I irecord, |
15 |
I myIter, |
16 |
I myThid ) |
17 |
C |
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 nNz integer size of third dimension: normally either 1 or Nr |
25 |
C arr RS/RL array to write, arr(:,:,nNz,:,:) |
26 |
C irecord integer record number to read |
27 |
C myIter integer time step number |
28 |
C myThid integer thread identifier |
29 |
C |
30 |
C MDS_WRITETILE creates either a file of the form "fName.data" and |
31 |
C "fName.meta" if the logical flag "globalFile" is set true. Otherwise |
32 |
C it creates MDS tiled files of the form "fName.xxx.yyy.data" and |
33 |
C "fName.xxx.yyy.meta". A meta-file is always created. |
34 |
C Currently, the meta-files are not read because it is difficult |
35 |
C to parse files in fortran. We should read meta information before |
36 |
C adding records to an existing multi-record file. |
37 |
C The precision of the file is decsribed by filePrec, set either |
38 |
C to floatPrec32 or floatPrec64. The precision or declaration of |
39 |
C the array argument must be consistently described by the char*(2) |
40 |
C string arrType, either "RS" or "RL". nNz allows for both 2-D and |
41 |
C 3-D arrays to be handled. nNz=1 implies a 2-D model field and |
42 |
C nNz=Nr implies a 3-D model field. irecord is the record number |
43 |
C to be read and must be >= 1. NOTE: It is currently assumed that |
44 |
C the highest record number in the file was the last record written. |
45 |
C Nor is there a consistency check between the routine arguments and file. |
46 |
C ie. if your write record 2 after record 4 the meta information |
47 |
C will record the number of records to be 2. This, again, is because |
48 |
C we have read the meta information. To be fixed. |
49 |
C |
50 |
C Created: 03/16/99 adcroft@mit.edu |
51 |
C |
52 |
C Changed: 05/31/00 heimbach@mit.edu |
53 |
C open(dUnit, ..., status='old', ... -> status='unknown' |
54 |
|
55 |
implicit none |
56 |
C Global variables / common blocks |
57 |
#include "SIZE.h" |
58 |
#include "EEPARAMS.h" |
59 |
#include "PARAMS.h" |
60 |
|
61 |
C Routine arguments |
62 |
character*(*) fName |
63 |
integer filePrec |
64 |
logical globalFile |
65 |
character*(2) arrType |
66 |
integer nNz |
67 |
_RL arr(1-oLx:sNx+oLx,1-oLy:sNy+oLy,nNz,nSx,nSy) |
68 |
integer bi, bj |
69 |
integer irecord |
70 |
integer myIter |
71 |
integer myThid |
72 |
C Functions |
73 |
integer ILNBLNK |
74 |
integer MDS_RECLEN |
75 |
C Local variables |
76 |
character*(80) dataFName,metaFName |
77 |
integer iG,jG,irec,j,k,dUnit,IL |
78 |
Real*4 r4seg(sNx) |
79 |
Real*8 r8seg(sNx) |
80 |
integer dimList(3,3),ndims |
81 |
integer length_of_rec |
82 |
logical fileIsOpen |
83 |
character*(max_len_mbuf) msgbuf |
84 |
C ------------------------------------------------------------------ |
85 |
|
86 |
C Only do I/O if I am the master thread |
87 |
_BEGIN_MASTER( myThid ) |
88 |
|
89 |
C Record number must be >= 1 |
90 |
if (irecord .LT. 1) then |
91 |
write(msgbuf,'(a,i9.8)') |
92 |
& ' MDS_WRITETILE: argument irecord = ',irecord |
93 |
call print_message( msgbuf, standardmessageunit, |
94 |
& SQUEEZE_RIGHT , mythid) |
95 |
write(msgbuf,'(a)') |
96 |
& ' MDS_WRITETILE: invalid value for irecord' |
97 |
call print_error( msgbuf, mythid ) |
98 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
99 |
endif |
100 |
|
101 |
C Assume nothing |
102 |
fileIsOpen=.FALSE. |
103 |
IL=ILNBLNK( fName ) |
104 |
|
105 |
C Assign a free unit number as the I/O channel for this routine |
106 |
call MDSFINDUNIT( dUnit, mythid ) |
107 |
|
108 |
C If we are writing to a global file then we open it here |
109 |
if (globalFile) then |
110 |
write(dataFname(1:80),'(2a)') fName(1:IL),'.data' |
111 |
if (irecord .EQ. 1) then |
112 |
length_of_rec=MDS_RECLEN( filePrec, sNx, mythid ) |
113 |
open( dUnit, file=dataFName, status='unknown', |
114 |
& access='direct', recl=length_of_rec ) |
115 |
fileIsOpen=.TRUE. |
116 |
else |
117 |
length_of_rec=MDS_RECLEN( filePrec, sNx, mythid ) |
118 |
open( dUnit, file=dataFName, status=_OLD_STATUS, |
119 |
& access='direct', recl=length_of_rec ) |
120 |
fileIsOpen=.TRUE. |
121 |
endif |
122 |
endif |
123 |
|
124 |
C Loop over all tiles |
125 |
c do bj=1,nSy |
126 |
c do bi=1,nSx |
127 |
C If we are writing to a tiled MDS file then we open each one here |
128 |
if (.NOT. globalFile) then |
129 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
130 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
131 |
write(dataFname(1:80),'(2a,i3.3,a,i3.3,a)') |
132 |
& fName(1:IL),'.',iG,'.',jG,'.data' |
133 |
if (irecord .EQ. 1) then |
134 |
length_of_rec=MDS_RECLEN( filePrec, sNx, mythid ) |
135 |
open( dUnit, file=dataFName, status=_NEW_STATUS, |
136 |
& access='direct', recl=length_of_rec ) |
137 |
fileIsOpen=.TRUE. |
138 |
else |
139 |
length_of_rec=MDS_RECLEN( filePrec, sNx, mythid ) |
140 |
open( dUnit, file=dataFName, status=_OLD_STATUS, |
141 |
& access='direct', recl=length_of_rec ) |
142 |
fileIsOpen=.TRUE. |
143 |
endif |
144 |
endif |
145 |
if (fileIsOpen) then |
146 |
do k=1,nNz |
147 |
do j=1,sNy |
148 |
if (globalFile) then |
149 |
iG = myXGlobalLo-1+(bi-1)*sNx |
150 |
jG = myYGlobalLo-1+(bj-1)*sNy |
151 |
irec=1+INT(iG/sNx)+nSx*nPx*(jG+j-1)+nSx*nPx*Ny*(k-1) |
152 |
& +nSx*nPx*Ny*nNz*(irecord-1) |
153 |
else |
154 |
iG = 0 |
155 |
jG = 0 |
156 |
irec=j + sNy*(k-1) + sNy*nNz*(irecord-1) |
157 |
endif |
158 |
if (filePrec .eq. precFloat32) then |
159 |
if (arrType .eq. 'RS') then |
160 |
call MDS_SEG4toRS( j,bi,bj,k,nNz, r4seg, .FALSE., arr ) |
161 |
elseif (arrType .eq. 'RL') then |
162 |
call MDS_SEG4toRL( j,bi,bj,k,nNz, r4seg, .FALSE., arr ) |
163 |
else |
164 |
write(msgbuf,'(a)') |
165 |
& ' MDS_WRITETILE: illegal value for arrType' |
166 |
call print_error( msgbuf, mythid ) |
167 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
168 |
endif |
169 |
#ifdef _BYTESWAPIO |
170 |
call MDS_BYTESWAPR4( sNx, r4seg ) |
171 |
#endif |
172 |
write(dUnit,rec=irec) r4seg |
173 |
elseif (filePrec .eq. precFloat64) then |
174 |
if (arrType .eq. 'RS') then |
175 |
call MDS_SEG8toRS( j,bi,bj,k,nNz, r8seg, .FALSE., arr ) |
176 |
elseif (arrType .eq. 'RL') then |
177 |
call MDS_SEG8toRL( j,bi,bj,k,nNz, r8seg, .FALSE., arr ) |
178 |
else |
179 |
write(msgbuf,'(a)') |
180 |
& ' MDS_WRITETILE: illegal value for arrType' |
181 |
call print_error( msgbuf, mythid ) |
182 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
183 |
endif |
184 |
#ifdef _BYTESWAPIO |
185 |
call MDS_BYTESWAPR8( sNx, r8seg ) |
186 |
#endif |
187 |
write(dUnit,rec=irec) r8seg |
188 |
else |
189 |
write(msgbuf,'(a)') |
190 |
& ' MDS_WRITETILE: illegal value for filePrec' |
191 |
call print_error( msgbuf, mythid ) |
192 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
193 |
endif |
194 |
C End of j loop |
195 |
enddo |
196 |
C End of k loop |
197 |
enddo |
198 |
else |
199 |
write(msgbuf,'(a)') |
200 |
& ' MDS_WRITETILE: I should never get to this point' |
201 |
call print_error( msgbuf, mythid ) |
202 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
203 |
endif |
204 |
C If we were writing to a tiled MDS file then we close it here |
205 |
if (fileIsOpen .AND. (.NOT. globalFile)) then |
206 |
close( dUnit ) |
207 |
fileIsOpen = .FALSE. |
208 |
endif |
209 |
C Create meta-file for each tile if we are tiling |
210 |
if (.NOT. globalFile) then |
211 |
iG=bi+(myXGlobalLo-1)/sNx ! Kludge until unstructered tiles |
212 |
jG=bj+(myYGlobalLo-1)/sNy ! Kludge until unstructered tiles |
213 |
write(metaFname(1:80),'(2a,i3.3,a,i3.3,a)') |
214 |
& fName(1:IL),'.',iG,'.',jG,'.meta' |
215 |
dimList(1,1)=Nx |
216 |
dimList(2,1)=myXGlobalLo+(bi-1)*sNx |
217 |
dimList(3,1)=myXGlobalLo+bi*sNx-1 |
218 |
dimList(1,2)=Ny |
219 |
dimList(2,2)=myYGlobalLo+(bj-1)*sNy |
220 |
dimList(3,2)=myYGlobalLo+bj*sNy-1 |
221 |
dimList(1,3)=Nr |
222 |
dimList(2,3)=1 |
223 |
dimList(3,3)=Nr |
224 |
ndims=3 |
225 |
if (nNz .EQ. 1) ndims=2 |
226 |
call MDSWRITEMETA( metaFName, dataFName, |
227 |
& filePrec, ndims, dimList, irecord, myIter, mythid ) |
228 |
endif |
229 |
C End of bi,bj loops |
230 |
c enddo |
231 |
c enddo |
232 |
|
233 |
C If global file was opened then close it |
234 |
if (fileIsOpen .AND. globalFile) then |
235 |
close( dUnit ) |
236 |
fileIsOpen = .FALSE. |
237 |
endif |
238 |
|
239 |
C Create meta-file for the global-file |
240 |
if (globalFile) then |
241 |
C We can not do this operation using threads (yet) because of the |
242 |
C "barrier" at the next step. The barrier could be removed but |
243 |
C at the cost of "safe" distributed I/O. |
244 |
if (nThreads.NE.1) then |
245 |
write(msgbuf,'(a,a)') |
246 |
& ' MDS_WRITETILE: A threads version of this routine', |
247 |
& ' does not exist.' |
248 |
call print_message( msgbuf, standardmessageunit, |
249 |
& SQUEEZE_RIGHT , mythid) |
250 |
write(msgbuf,'(a)') |
251 |
& ' MDS_WRITETILE: This needs to be fixed...' |
252 |
call print_message( msgbuf, standardmessageunit, |
253 |
& SQUEEZE_RIGHT , mythid) |
254 |
write(msgbuf,'(a,i3.2)') |
255 |
& ' MDS_WRITETILE: nThreads = ',nThreads |
256 |
call print_message( msgbuf, standardmessageunit, |
257 |
& SQUEEZE_RIGHT , mythid) |
258 |
write(msgbuf,'(a)') |
259 |
& ' MDS_WRITETILE: Stopping because you are using threads' |
260 |
call print_error( msgbuf, mythid ) |
261 |
stop 'ABNORMAL END: S/R MDS_WRITETILE' |
262 |
endif |
263 |
C We put a barrier here to ensure that all processes have finished |
264 |
C writing their data before we update the meta-file |
265 |
_BARRIER |
266 |
write(metaFName(1:80),'(2a)') fName(1:IL),'.meta' |
267 |
dimList(1,1)=Nx |
268 |
dimList(2,1)=1 |
269 |
dimList(3,1)=Nx |
270 |
dimList(1,2)=Ny |
271 |
dimList(2,2)=1 |
272 |
dimList(3,2)=Ny |
273 |
dimList(1,3)=Nr |
274 |
dimList(2,3)=1 |
275 |
dimList(3,3)=Nr |
276 |
ndims=3 |
277 |
if (nNz .EQ. 1) ndims=2 |
278 |
call MDSWRITEMETA( metaFName, dataFName, |
279 |
& filePrec, ndims, dimList, irecord, myIter, mythid ) |
280 |
fileIsOpen=.TRUE. |
281 |
endif |
282 |
|
283 |
_END_MASTER( myThid ) |
284 |
|
285 |
C ------------------------------------------------------------------ |
286 |
return |
287 |
end |