1 |
C $Header: /u/gcmpack/MITgcm/pkg/mnc/mnc_cw_readwrite.template,v 1.22 2004/10/20 21:26:14 edhill Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "MNC_OPTIONS.h" |
5 |
|
6 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
7 |
CBOP 0 |
8 |
C !ROUTINE: MNC_CW_RX_W_S |
9 |
|
10 |
C !INTERFACE: |
11 |
SUBROUTINE MNC_CW_RX_W_S( |
12 |
I stype, |
13 |
I fbname, bi,bj, |
14 |
I vtype, |
15 |
I var, |
16 |
I myThid ) |
17 |
|
18 |
C !DESCRIPTION: |
19 |
C A scalar version of MNC_CW_RX_W() for compilers that cannot |
20 |
C gracefully handle the conversion on their own. |
21 |
|
22 |
C !USES: |
23 |
implicit none |
24 |
|
25 |
C !INPUT PARAMETERS: |
26 |
integer myThid, bi,bj |
27 |
character*(*) stype, fbname, vtype |
28 |
__V var |
29 |
__V var_arr(1) |
30 |
CEOP |
31 |
|
32 |
var_arr(1) = var |
33 |
CALL MNC_CW_RX_W(stype,fbname,bi,bj,vtype, var_arr, myThid) |
34 |
|
35 |
RETURN |
36 |
END |
37 |
|
38 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
39 |
CBOP 0 |
40 |
C !ROUTINE: MNC_CW_RX_W |
41 |
|
42 |
C !INTERFACE: |
43 |
SUBROUTINE MNC_CW_RX_W( |
44 |
I stype, |
45 |
I fbname, bi,bj, |
46 |
I vtype, |
47 |
I var, |
48 |
I myThid ) |
49 |
|
50 |
C !DESCRIPTION: |
51 |
C A scalar version of MNC_CW_RX_W() for compilers that cannot |
52 |
C gracefully handle the conversion on their own. |
53 |
|
54 |
C !USES: |
55 |
implicit none |
56 |
|
57 |
C !INPUT PARAMETERS: |
58 |
integer myThid, bi,bj |
59 |
character*(*) stype, fbname, vtype |
60 |
__V var(*) |
61 |
INTEGER offsets(9) |
62 |
CEOP |
63 |
INTEGER i |
64 |
|
65 |
DO i = 1,9 |
66 |
offsets(i) = 0 |
67 |
ENDDO |
68 |
CALL MNC_CW_RX_W_OFFSET(stype,fbname,bi,bj,vtype, var, |
69 |
& offsets, myThid) |
70 |
|
71 |
RETURN |
72 |
END |
73 |
|
74 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
75 |
CBOP 0 |
76 |
C !ROUTINE: MNC_CW_RX_W_OFFSET |
77 |
|
78 |
C !INTERFACE: |
79 |
SUBROUTINE MNC_CW_RX_W_OFFSET( |
80 |
I stype, |
81 |
I fbname, bi,bj, |
82 |
I vtype, |
83 |
I var, |
84 |
I offsets, |
85 |
I myThid ) |
86 |
|
87 |
C !DESCRIPTION: |
88 |
C This subroutine writes one variable to a file or a file group, |
89 |
C depending upon the tile indicies. |
90 |
|
91 |
C !USES: |
92 |
implicit none |
93 |
#include "netcdf.inc" |
94 |
#include "mnc_common.h" |
95 |
#include "SIZE.h" |
96 |
#include "EEPARAMS.h" |
97 |
#include "PARAMS.h" |
98 |
#include "MNC_PARAMS.h" |
99 |
|
100 |
C !INPUT PARAMETERS: |
101 |
integer myThid, bi,bj |
102 |
character*(*) stype, fbname, vtype |
103 |
__V var(*) |
104 |
INTEGER offsets(*) |
105 |
CEOP |
106 |
|
107 |
C !LOCAL VARIABLES: |
108 |
integer i,j,k, indv,nvf,nvl, n1,n2, igrid, ntot, indu |
109 |
integer bis,bie, bjs,bje, uniq_tnum, nfname, iseq |
110 |
integer fid, idv, indvids, ndim, indf, err |
111 |
integer lbi,lbj, bidim,bjdim, unlim_sz, kr |
112 |
integer p(9),s(9),e(9), dimnc(9) |
113 |
integer vstart(9),vcount(9), udo(9) |
114 |
integer j1,j2,j3,j4,j5,j6,j7, k1,k2,k3,k4,k5,k6,k7 |
115 |
integer indfg, fg1,fg2, npath |
116 |
character*(MAX_LEN_MBUF) msgbuf |
117 |
character*(MNC_MAX_CHAR) fname |
118 |
character*(MNC_MAX_CHAR) path_fname |
119 |
REAL*8 resh_d( sNx + 2*OLx + sNy + 2*OLy ) |
120 |
REAL*4 resh_r( sNx + 2*OLx + sNy + 2*OLy ) |
121 |
INTEGER resh_i( sNx + 2*OLx + sNy + 2*OLy ) |
122 |
#ifdef HAVE_STAT |
123 |
integer ntotenc, ncenc, nbytes, fs_isdone |
124 |
character*(200) cenc |
125 |
integer ienc(200) |
126 |
REAL*8 fsnu |
127 |
#endif |
128 |
|
129 |
C Functions |
130 |
integer IFNBLNK, ILNBLNK |
131 |
|
132 |
C Only do I/O if I am the master thread |
133 |
_BEGIN_MASTER( myThid ) |
134 |
|
135 |
C Get the current index for the unlimited dimension from the file |
136 |
C group (or base) name |
137 |
fg1 = IFNBLNK(fbname) |
138 |
fg2 = ILNBLNK(fbname) |
139 |
CALL MNC_GET_IND(MNC_MAX_ID, fbname, mnc_cw_fgnm, indfg, myThid) |
140 |
IF (indfg .LT. 1) THEN |
141 |
write(msgbuf,'(3a)') |
142 |
& 'MNC_CW_RX_W ERROR: file group name ''', |
143 |
& fbname(fg1:fg2), ''' is not defined' |
144 |
CALL print_error(msgbuf, mythid) |
145 |
STOP 'ABNORMAL END: S/R MNC_CW_RX_W' |
146 |
ENDIF |
147 |
indu = mnc_cw_fgud(indfg) |
148 |
iseq = mnc_cw_fgis(indfg) |
149 |
C write(*,*) 'indu,iseq = ', indu, iseq |
150 |
|
151 |
C Check that the Variable Type exists |
152 |
nvf = IFNBLNK(vtype) |
153 |
nvl = ILNBLNK(vtype) |
154 |
CALL MNC_GET_IND(MNC_MAX_ID, vtype, mnc_cw_vname, indv, myThid) |
155 |
IF (indv .LT. 1) THEN |
156 |
write(msgbuf,'(3a)') 'MNC_CW_RX_W ERROR: vtype ''', |
157 |
& vtype(nvf:nvl), ''' is not defined' |
158 |
CALL print_error(msgbuf, mythid) |
159 |
STOP 'ABNORMAL END: S/R MNC_CW_RX_W' |
160 |
ENDIF |
161 |
igrid = mnc_cw_vgind(indv) |
162 |
|
163 |
C Set the bi,bj indicies |
164 |
bis = bi |
165 |
bie = bi |
166 |
IF (bi .LT. 1) THEN |
167 |
bis = 1 |
168 |
bie = nSx |
169 |
ENDIF |
170 |
bjs = bj |
171 |
bje = bj |
172 |
IF (bj .LT. 1) THEN |
173 |
bjs = 1 |
174 |
bje = nSy |
175 |
ENDIF |
176 |
|
177 |
DO lbj = bjs,bje |
178 |
DO lbi = bis,bie |
179 |
|
180 |
#ifdef HAVE_STAT |
181 |
fs_isdone = 0 |
182 |
#endif |
183 |
10 CONTINUE |
184 |
|
185 |
C Create the file name |
186 |
CALL MNC_CW_GET_TILE_NUM(lbi,lbj, uniq_tnum, myThid) |
187 |
fname(1:MNC_MAX_CHAR) = mnc_blank_name(1:MNC_MAX_CHAR) |
188 |
n1 = IFNBLNK(fbname) |
189 |
n2 = ILNBLNK(fbname) |
190 |
ntot = n2 - n1 + 1 |
191 |
fname(1:ntot) = fbname(n1:n2) |
192 |
ntot = ntot + 1 |
193 |
fname(ntot:ntot) = '.' |
194 |
write(fname((ntot+1):(ntot+14)),'(i4.4,a1,i6.6,a3)') |
195 |
& iseq,'.',uniq_tnum, '.nc' |
196 |
nfname = ntot + 14 |
197 |
|
198 |
C Add the path to the file name |
199 |
IF (mnc_use_outdir) THEN |
200 |
path_fname(1:MNC_MAX_CHAR) = mnc_blank_name(1:MNC_MAX_CHAR) |
201 |
npath = ILNBLNK(mnc_out_path) |
202 |
path_fname(1:npath) = mnc_out_path(1:npath) |
203 |
path_fname((npath+1):(npath+nfname)) = fname(1:nfname) |
204 |
fname(1:MNC_MAX_CHAR) = path_fname(1:MNC_MAX_CHAR) |
205 |
nfname = npath + nfname |
206 |
ENDIF |
207 |
|
208 |
C Append to an existing or create a new file |
209 |
CALL MNC_CW_FILE_AORC(fname, indf, myThid) |
210 |
fid = mnc_f_info(indf,2) |
211 |
|
212 |
#ifdef HAVE_STAT |
213 |
IF ((mnc_cw_fgig(indfg) .EQ. 1) |
214 |
& .AND. (fs_isdone .EQ. 0)) THEN |
215 |
C Decide whether to append to the existing or create a new file |
216 |
C based on the byte count per unlimited dimension |
217 |
ncenc = 70 |
218 |
cenc(1:26) = 'abcdefghijklmnopqrstuvwxyz' |
219 |
cenc(27:52) = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ' |
220 |
cenc(53:70) = '0123456789_.,+-=/~' |
221 |
k = nfname |
222 |
IF (k .GT. 200) k = 200 |
223 |
ntotenc = 0 |
224 |
DO i = 1,k |
225 |
DO j = 1,ncenc |
226 |
IF (fname(i:i) .EQ. cenc(j:j)) THEN |
227 |
ntotenc = ntotenc + 1 |
228 |
ienc(ntotenc) = j |
229 |
GOTO 20 |
230 |
ENDIF |
231 |
ENDDO |
232 |
20 CONTINUE |
233 |
ENDDO |
234 |
CALL mncfsize(ntotenc, ienc, nbytes) |
235 |
IF (nbytes .GT. 0) THEN |
236 |
CALL MNC_DIM_UNLIM_SIZE(fname, unlim_sz, myThid) |
237 |
fsnu = (1.0 _d 0 + 1.0 _d 0 / DBLE(unlim_sz)) |
238 |
& * DBLE(nbytes) |
239 |
IF (fsnu .GT. mnc_max_fsize) THEN |
240 |
iseq = iseq + 1 |
241 |
indu = 1 |
242 |
mnc_cw_fgud(indfg) = 1 |
243 |
mnc_cw_fgis(indfg) = iseq |
244 |
fs_isdone = 1 |
245 |
GOTO 10 |
246 |
ENDIF |
247 |
ENDIF |
248 |
ENDIF |
249 |
#endif /* HAVE_STAT */ |
250 |
|
251 |
C Ensure that all the NetCDF dimensions are defined and create a |
252 |
C local copy of them |
253 |
DO i = 1,9 |
254 |
dimnc(i) = 1 |
255 |
ENDDO |
256 |
DO i = 1,mnc_cw_ndim(igrid) |
257 |
IF (mnc_cw_dims(i,igrid) .EQ. -1) THEN |
258 |
dimnc(i) = -1 |
259 |
ELSE |
260 |
dimnc(i) = mnc_cw_ie(i,igrid) - mnc_cw_is(i,igrid) + 1 |
261 |
ENDIF |
262 |
CALL MNC_DIM_INIT(fname, |
263 |
& mnc_cw_dn(i,igrid), dimnc(i), myThid) |
264 |
ENDDO |
265 |
|
266 |
C Ensure that the "grid" is defined |
267 |
CALL MNC_GRID_INIT(fname, mnc_cw_gname(igrid), |
268 |
& mnc_cw_ndim(igrid), mnc_cw_dn(1,igrid), myThid) |
269 |
|
270 |
C Ensure that the variable is defined |
271 |
IF (stype(1:1) .EQ. 'D') |
272 |
& CALL MNC_VAR_INIT_DBL( |
273 |
& fname, mnc_cw_gname(igrid), vtype, myThid) |
274 |
IF (stype(1:1) .EQ. 'R') |
275 |
& CALL MNC_VAR_INIT_REAL( |
276 |
& fname, mnc_cw_gname(igrid), vtype, myThid) |
277 |
IF (stype(1:1) .EQ. 'I') |
278 |
& CALL MNC_VAR_INIT_INT( |
279 |
& fname, mnc_cw_gname(igrid), vtype, myThid) |
280 |
|
281 |
DO i = 1,mnc_fv_ids(indf,1) |
282 |
j = 2 + 3*(i - 1) |
283 |
IF (mnc_v_names(mnc_fv_ids(indf,j)) .EQ. vtype) THEN |
284 |
idv = mnc_fv_ids(indf,j+1) |
285 |
indvids = mnc_fd_ind(indf, mnc_f_info(indf, |
286 |
& (mnc_fv_ids(indf,j+2) + 1)) ) |
287 |
GOTO 30 |
288 |
ENDIF |
289 |
ENDDO |
290 |
write(msgbuf,'(4a)') 'MNC_MNC_CW_RX_W ERROR: ', |
291 |
& 'cannot reference variable ''', vtype, '''' |
292 |
CALL print_error(msgbuf, mythid) |
293 |
STOP 'ABNORMAL END: package MNC' |
294 |
30 CONTINUE |
295 |
|
296 |
C Check for bi,bj indicies |
297 |
bidim = mnc_cw_vbij(1,indv) |
298 |
bjdim = mnc_cw_vbij(2,indv) |
299 |
CEH3 write(*,*) 'bidim,bjdim = ', bidim,bjdim |
300 |
|
301 |
C Set the dimensions for the in-memory array |
302 |
ndim = mnc_cw_ndim(igrid) |
303 |
k = mnc_cw_dims(1,igrid) |
304 |
IF (k .GT. 0) THEN |
305 |
p(1) = k |
306 |
ELSE |
307 |
p(1) = 1 |
308 |
ENDIF |
309 |
DO i = 2,9 |
310 |
k = mnc_cw_dims(i,igrid) |
311 |
IF (k .LT. 1) THEN |
312 |
k = 1 |
313 |
ENDIF |
314 |
IF ((bidim .GT. 0) .AND. (i .EQ. bidim)) THEN |
315 |
p(i) = nSx * p(i-1) |
316 |
ELSEIF ((bjdim .GT. 0) .AND. (i .EQ. bjdim)) THEN |
317 |
p(i) = nSy * p(i-1) |
318 |
ELSE |
319 |
p(i) = k * p(i-1) |
320 |
ENDIF |
321 |
IF (offsets(i) .GT. 0) THEN |
322 |
k = 1 |
323 |
p(i) = k * p(i-1) |
324 |
ENDIF |
325 |
ENDDO |
326 |
|
327 |
C Set starting and ending indicies for the in-memory array and |
328 |
C the unlimited dimension offset for the NetCDF array |
329 |
DO i = 1,9 |
330 |
udo(i) = 0 |
331 |
s(i) = 1 |
332 |
e(i) = 1 |
333 |
IF (i .LE. ndim) THEN |
334 |
s(i) = mnc_cw_is(i,igrid) |
335 |
e(i) = mnc_cw_ie(i,igrid) |
336 |
ENDIF |
337 |
C Check for the unlimited dimension |
338 |
IF ((i .EQ. ndim) |
339 |
& .AND. (mnc_cw_dims(i,igrid) .EQ. -1)) THEN |
340 |
IF (indu .GT. 0) THEN |
341 |
C Use the indu value |
342 |
udo(i) = indu - 1 |
343 |
ELSEIF (indu .EQ. -1) THEN |
344 |
C Append one to the current unlimited dim size |
345 |
CALL MNC_DIM_UNLIM_SIZE( fname, unlim_sz, myThid) |
346 |
udo(i) = unlim_sz |
347 |
ELSE |
348 |
C Use the current unlimited dim size |
349 |
CALL MNC_DIM_UNLIM_SIZE( fname, unlim_sz, myThid) |
350 |
udo(i) = unlim_sz - 1 |
351 |
ENDIF |
352 |
ENDIF |
353 |
ENDDO |
354 |
IF (bidim .GT. 0) THEN |
355 |
s(bidim) = lbi |
356 |
e(bidim) = lbi |
357 |
ENDIF |
358 |
IF (bjdim .GT. 0) THEN |
359 |
s(bjdim) = lbj |
360 |
e(bjdim) = lbj |
361 |
ENDIF |
362 |
|
363 |
C Check the offsets |
364 |
DO i = 1,9 |
365 |
IF (offsets(i) .GT. 0) THEN |
366 |
udo(i) = udo(i) + offsets(i) - 1 |
367 |
s(i) = 1 |
368 |
e(i) = 1 |
369 |
ENDIF |
370 |
ENDDO |
371 |
|
372 |
C Add the global attributes |
373 |
CALL MNC_CW_SET_GATTR( fname, lbi,lbj, uniq_tnum, myThid) |
374 |
|
375 |
C Add the per-variable attributes |
376 |
DO i = 1,mnc_cw_vnat(1,indv) |
377 |
CALL MNC_VAR_ADD_ATTR_STR( fname, vtype, |
378 |
& mnc_cw_vtnm(i,indv), mnc_cw_vtat(i,indv), myThid) |
379 |
ENDDO |
380 |
DO i = 1,mnc_cw_vnat(2,indv) |
381 |
CALL MNC_VAR_ADD_ATTR_INT( fname, vtype, |
382 |
& mnc_cw_vinm(i,indv), 1, mnc_cw_viat(i,indv), myThid) |
383 |
ENDDO |
384 |
DO i = 1,mnc_cw_vnat(3,indv) |
385 |
CALL MNC_VAR_ADD_ATTR_DBL( fname, vtype, |
386 |
& mnc_cw_vdnm(i,indv), 1, mnc_cw_vdat(i,indv), myThid) |
387 |
ENDDO |
388 |
|
389 |
CALL MNC_FILE_ENDDEF(fname, myThid) |
390 |
|
391 |
write(msgbuf,'(5a)') 'writing variable type ''', |
392 |
& vtype(nvf:nvl), ''' within file ''', |
393 |
& fname(1:nfname), '''' |
394 |
|
395 |
C DO i = 1,9 |
396 |
C write(*,*) 'i,p(i),s(i),e(i),udo(i),offsets(i) = ', |
397 |
C & i,p(i),s(i),e(i),udo(i),offsets(i) |
398 |
C ENDDO |
399 |
|
400 |
C Write the variable one vector at a time |
401 |
DO j7 = s(7),e(7) |
402 |
k7 = (j7 - 1)*p(6) |
403 |
vstart(7) = udo(7) + j7 - s(7) + 1 |
404 |
vcount(7) = 1 |
405 |
DO j6 = s(6),e(6) |
406 |
k6 = (j6 - 1)*p(5) + k7 |
407 |
vstart(6) = udo(6) + j6 - s(6) + 1 |
408 |
vcount(6) = 1 |
409 |
DO j5 = s(5),e(5) |
410 |
k5 = (j5 - 1)*p(4) + k6 |
411 |
vstart(5) = udo(5) + j5 - s(5) + 1 |
412 |
vcount(5) = 1 |
413 |
DO j4 = s(4),e(4) |
414 |
k4 = (j4 - 1)*p(3) + k5 |
415 |
vstart(4) = udo(4) + j4 - s(4) + 1 |
416 |
vcount(4) = 1 |
417 |
DO j3 = s(3),e(3) |
418 |
k3 = (j3 - 1)*p(2) + k4 |
419 |
vstart(3) = udo(3) + j3 - s(3) + 1 |
420 |
vcount(3) = 1 |
421 |
DO j2 = s(2),e(2) |
422 |
k2 = (j2 - 1)*p(1) + k3 |
423 |
vstart(2) = udo(2) + j2 - s(2) + 1 |
424 |
vcount(2) = 1 |
425 |
|
426 |
kr = 0 |
427 |
vstart(1) = udo(1) + 1 |
428 |
vcount(1) = e(1) - s(1) + 1 |
429 |
|
430 |
IF (stype(1:1) .EQ. 'D') THEN |
431 |
DO j1 = s(1),e(1) |
432 |
k1 = k2 + j1 |
433 |
kr = kr + 1 |
434 |
resh_d(kr) = var(k1) |
435 |
ENDDO |
436 |
err = NF_PUT_VARA_DOUBLE(fid, idv, vstart, vcount, resh_d) |
437 |
ENDIF |
438 |
IF (stype(1:1) .EQ. 'R') THEN |
439 |
DO j1 = s(1),e(1) |
440 |
k1 = k2 + j1 |
441 |
kr = kr + 1 |
442 |
resh_r(kr) = var(k1) |
443 |
ENDDO |
444 |
err = NF_PUT_VARA_REAL(fid, idv, vstart, vcount, resh_r) |
445 |
ENDIF |
446 |
IF (stype(1:1) .EQ. 'I') THEN |
447 |
DO j1 = s(1),e(1) |
448 |
k1 = k2 + j1 |
449 |
kr = kr + 1 |
450 |
resh_i(kr) = var(k1) |
451 |
ENDDO |
452 |
err = NF_PUT_VARA_INT(fid, idv, vstart, vcount, resh_i) |
453 |
ENDIF |
454 |
|
455 |
CALL MNC_HANDLE_ERR(err, msgbuf, myThid) |
456 |
|
457 |
ENDDO |
458 |
ENDDO |
459 |
ENDDO |
460 |
ENDDO |
461 |
ENDDO |
462 |
ENDDO |
463 |
|
464 |
C Sync the file |
465 |
err = NF_SYNC(fid) |
466 |
write(msgbuf,'(3a)') 'sync for file ''', fname, |
467 |
& ''' in S/R MNC_CW_RX_W' |
468 |
CALL MNC_HANDLE_ERR(err, msgbuf, myThid) |
469 |
|
470 |
ENDDO |
471 |
ENDDO |
472 |
|
473 |
_END_MASTER( myThid ) |
474 |
|
475 |
RETURN |
476 |
END |
477 |
|
478 |
|
479 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
480 |
CBOP 0 |
481 |
C !ROUTINE: MNC_CW_RX_R_S |
482 |
|
483 |
C !INTERFACE: |
484 |
SUBROUTINE MNC_CW_RX_R_S( |
485 |
I stype, |
486 |
I fbname, bi,bj, |
487 |
I vtype, |
488 |
I var, |
489 |
I myThid ) |
490 |
|
491 |
C !DESCRIPTION: |
492 |
C A scalar version of MNC_CW_RX_R() for compilers that cannot |
493 |
C gracefully handle the conversion on their own. |
494 |
|
495 |
C !USES: |
496 |
implicit none |
497 |
|
498 |
C !INPUT PARAMETERS: |
499 |
integer myThid, bi,bj |
500 |
character*(*) stype, fbname, vtype |
501 |
__V var |
502 |
__V var_arr(1) |
503 |
CEOP |
504 |
var_arr(1) = var |
505 |
|
506 |
CALL MNC_CW_RX_R(stype,fbname,bi,bj,vtype, var_arr, myThid) |
507 |
|
508 |
RETURN |
509 |
END |
510 |
|
511 |
|
512 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
513 |
CBOP 0 |
514 |
C !ROUTINE: MNC_CW_RX_R |
515 |
|
516 |
C !INTERFACE: |
517 |
SUBROUTINE MNC_CW_RX_R( |
518 |
I stype, |
519 |
I fbname, bi,bj, |
520 |
I vtype, |
521 |
I var, |
522 |
I myThid ) |
523 |
|
524 |
C !DESCRIPTION: |
525 |
C This subroutine reads one variable from a file or a file group, |
526 |
C depending upon the tile indicies. |
527 |
|
528 |
C !USES: |
529 |
implicit none |
530 |
#include "netcdf.inc" |
531 |
#include "mnc_common.h" |
532 |
#include "SIZE.h" |
533 |
#include "EEPARAMS.h" |
534 |
#include "PARAMS.h" |
535 |
#include "MNC_PARAMS.h" |
536 |
|
537 |
C !INPUT PARAMETERS: |
538 |
integer myThid, bi,bj |
539 |
character*(*) stype, fbname, vtype |
540 |
__V var(*) |
541 |
CEOP |
542 |
|
543 |
C !LOCAL VARIABLES: |
544 |
integer i,k, nvf,nvl, n1,n2, igrid, ntot, indu |
545 |
integer bis,bie, bjs,bje, uniq_tnum, nfname, fid, idv |
546 |
integer ndim, indf, err, lbi,lbj, bidim,bjdim, unlim_sz, kr |
547 |
integer ind_fv_ids, ind_vt, ierr, atype, alen |
548 |
integer f_sNx,f_sNy, npath |
549 |
integer p(9),s(9),e(9), vstart(9),vcount(9), udo(9) |
550 |
integer j1,j2,j3,j4,j5,j6,j7, k1,k2,k3,k4,k5,k6,k7 |
551 |
character*(MAX_LEN_MBUF) msgbuf |
552 |
character*(MNC_MAX_CHAR) fname |
553 |
character*(MNC_MAX_CHAR) path_fname |
554 |
integer indfg, fg1,fg2 |
555 |
REAL*8 resh_d( sNx + 2*OLx + sNy + 2*OLy ) |
556 |
REAL*4 resh_r( sNx + 2*OLx + sNy + 2*OLy ) |
557 |
INTEGER resh_i( sNx + 2*OLx + sNy + 2*OLy ) |
558 |
|
559 |
C Functions |
560 |
integer IFNBLNK, ILNBLNK |
561 |
|
562 |
C Only do I/O if I am the master thread |
563 |
_BEGIN_MASTER( myThid ) |
564 |
|
565 |
C Get the current index for the unlimited dimension from the file |
566 |
C group (or base) name |
567 |
fg1 = IFNBLNK(fbname) |
568 |
fg2 = ILNBLNK(fbname) |
569 |
CALL MNC_GET_IND(MNC_MAX_ID, fbname, mnc_cw_fgnm, indfg, myThid) |
570 |
IF (indfg .LT. 1) THEN |
571 |
write(msgbuf,'(3a)') |
572 |
& 'MNC_CW_RX_W ERROR: file group name ''', |
573 |
& fbname(fg1:fg2), ''' is not defined' |
574 |
CALL print_error(msgbuf, mythid) |
575 |
STOP 'ABNORMAL END: S/R MNC_CW_RX_W' |
576 |
ENDIF |
577 |
indu = mnc_cw_fgud(indfg) |
578 |
|
579 |
C Check that the Variable Type exists |
580 |
nvf = IFNBLNK(vtype) |
581 |
nvl = ILNBLNK(vtype) |
582 |
CALL MNC_GET_IND( MNC_MAX_ID, vtype, mnc_cw_vname, ind_vt, myThid) |
583 |
IF (ind_vt .LT. 1) THEN |
584 |
write(msgbuf,'(3a)') 'MNC_CW_RX_R ERROR: vtype ''', |
585 |
& vtype(nvf:nvl), ''' is not defined' |
586 |
CALL print_error(msgbuf, mythid) |
587 |
STOP 'ABNORMAL END: S/R MNC_CW_RX_R' |
588 |
ENDIF |
589 |
igrid = mnc_cw_vgind(ind_vt) |
590 |
|
591 |
C Check for bi,bj indicies |
592 |
bidim = mnc_cw_vbij(1,ind_vt) |
593 |
bjdim = mnc_cw_vbij(2,ind_vt) |
594 |
|
595 |
C Set the bi,bj indicies |
596 |
bis = bi |
597 |
bie = bi |
598 |
IF (bi .LT. 1) THEN |
599 |
bis = 1 |
600 |
bie = nSx |
601 |
ENDIF |
602 |
bjs = bj |
603 |
bje = bj |
604 |
IF (bj .LT. 1) THEN |
605 |
bjs = 1 |
606 |
bje = nSy |
607 |
ENDIF |
608 |
|
609 |
DO lbj = bjs,bje |
610 |
DO lbi = bis,bie |
611 |
|
612 |
C Create the file name |
613 |
CALL MNC_CW_GET_TILE_NUM( lbi,lbj, uniq_tnum, myThid) |
614 |
fname(1:MNC_MAX_CHAR) = mnc_blank_name(1:MNC_MAX_CHAR) |
615 |
n1 = IFNBLNK(fbname) |
616 |
n2 = ILNBLNK(fbname) |
617 |
ntot = n2 - n1 + 1 |
618 |
fname(1:ntot) = fbname(n1:n2) |
619 |
ntot = ntot + 1 |
620 |
fname(ntot:ntot) = '.' |
621 |
write(fname((ntot+1):(ntot+9)),'(i6.6,a3)') uniq_tnum, '.nc' |
622 |
nfname = ntot+9 |
623 |
|
624 |
C Add the path to the file name |
625 |
IF (mnc_use_indir) THEN |
626 |
path_fname(1:MNC_MAX_CHAR) = mnc_blank_name(1:MNC_MAX_CHAR) |
627 |
npath = ILNBLNK(mnc_indir_str) |
628 |
path_fname(1:npath) = mnc_indir_str(1:npath) |
629 |
path_fname((npath+1):(npath+nfname)) = fname(1:nfname) |
630 |
fname(1:MNC_MAX_CHAR) = path_fname(1:MNC_MAX_CHAR) |
631 |
nfname = npath + nfname |
632 |
ENDIF |
633 |
|
634 |
C Open the existing file |
635 |
CALL MNC_FILE_TRY_READ( fname, ierr, indf, myThid) |
636 |
|
637 |
C Check that the variable (VType) is defined within the file |
638 |
CALL MNC_GET_FVINDS( fname, vtype, indf, ind_fv_ids, myThid) |
639 |
IF ((indf .LT. 1) .OR. (ind_fv_ids .LT. 1)) THEN |
640 |
write(msgbuf,'(4a)') 'MNC_CW_RX_R ERROR: vtype ''', |
641 |
& vtype(nvf:nvl), ''' is not defined within file ''', |
642 |
& fname(1:nfname) |
643 |
CALL print_error(msgbuf, mythid) |
644 |
STOP 'ABNORMAL END: S/R MNC_CW_RX_R' |
645 |
ENDIF |
646 |
fid = mnc_f_info(indf,2) |
647 |
idv = mnc_fv_ids(indf,ind_fv_ids+1) |
648 |
|
649 |
C Check that the current sNy,sNy values and the in-file values |
650 |
C are compatible and WARN (only warn) if not |
651 |
f_sNx = -1 |
652 |
f_sNy = -1 |
653 |
err = NF_INQ_ATT(fid,NF_GLOBAL, 'sNx',atype,alen) |
654 |
IF ((err .EQ. NF_NOERR) .AND. (alen .EQ. 1)) THEN |
655 |
err = NF_GET_ATT_INT(fid, NF_GLOBAL, 'sNx', f_sNx) |
656 |
CALL MNC_HANDLE_ERR(err, |
657 |
& 'reading attribute ''sNx'' in S/R MNC_CW_RX_R', |
658 |
& myThid) |
659 |
ENDIF |
660 |
err = NF_INQ_ATT(fid,NF_GLOBAL, 'sNy',atype,alen) |
661 |
IF ((err .EQ. NF_NOERR) .AND. (alen .EQ. 1)) THEN |
662 |
err = NF_GET_ATT_INT(fid, NF_GLOBAL, 'sNy', f_sNy) |
663 |
CALL MNC_HANDLE_ERR(err, |
664 |
& 'reading attribute ''sNy'' in S/R MNC_CW_RX_R', |
665 |
& myThid) |
666 |
ENDIF |
667 |
IF ((f_sNx .NE. sNx) .OR. (f_sNy .NE. sNy)) THEN |
668 |
write(msgbuf,'(5a)') 'MNC_CW_RX_R WARNING: the ', |
669 |
& 'attributes ''sNx'' and ''sNy'' within the file ''', |
670 |
& fname(1:nfname), ''' do not exist or do not match ', |
671 |
& 'the current sizes within the model' |
672 |
CALL print_error(msgbuf, mythid) |
673 |
ENDIF |
674 |
|
675 |
C Check that the in-memory variable and the in-file variables |
676 |
C are of compatible sizes |
677 |
C ires = 1 |
678 |
C CALL MNC_CHK_VTYP_R_NCVAR( ind_vt, |
679 |
C & indf, ind_fv_ids, indu, ires) |
680 |
C IF (ires .LT. 0) THEN |
681 |
C write(msgbuf,'(7a)') 'MNC_CW_RX_R WARNING: the sizes ', |
682 |
C & 'of the in-program variable ''', vtype(nvf:nvl), |
683 |
C & ''' and the corresponding variable within file ''', |
684 |
C & fname(1:nfname), ''' are not compatible -- please ', |
685 |
C & 'check the sizes' |
686 |
C CALL print_error(msgbuf, mythid) |
687 |
C STOP 'ABNORMAL END: S/R MNC_CW_RX_R' |
688 |
C ENDIF |
689 |
|
690 |
C Check for bi,bj indicies |
691 |
bidim = mnc_cw_vbij(1,ind_vt) |
692 |
bjdim = mnc_cw_vbij(2,ind_vt) |
693 |
|
694 |
C Set the dimensions for the in-memory array |
695 |
ndim = mnc_cw_ndim(igrid) |
696 |
k = mnc_cw_dims(1,igrid) |
697 |
IF (k .GT. 0) THEN |
698 |
p(1) = k |
699 |
ELSE |
700 |
p(1) = 1 |
701 |
ENDIF |
702 |
DO i = 2,9 |
703 |
k = mnc_cw_dims(i,igrid) |
704 |
IF (k .LT. 1) THEN |
705 |
k = 1 |
706 |
ENDIF |
707 |
IF ((bidim .GT. 0) .AND. (i .EQ. bidim)) THEN |
708 |
p(i) = nSx * p(i-1) |
709 |
ELSEIF ((bjdim .GT. 0) .AND. (i .EQ. bjdim)) THEN |
710 |
p(i) = nSy * p(i-1) |
711 |
ELSE |
712 |
p(i) = k * p(i-1) |
713 |
ENDIF |
714 |
ENDDO |
715 |
|
716 |
C Set starting and ending indicies for the in-memory array and |
717 |
C the unlimited dimension offset for the NetCDF array |
718 |
DO i = 1,9 |
719 |
udo(i) = 0 |
720 |
s(i) = 1 |
721 |
e(i) = 1 |
722 |
IF (i .LE. ndim) THEN |
723 |
s(i) = mnc_cw_is(i,igrid) |
724 |
e(i) = mnc_cw_ie(i,igrid) |
725 |
ENDIF |
726 |
C Check for the unlimited dimension |
727 |
IF ((i .EQ. ndim) |
728 |
& .AND. (mnc_cw_dims(i,igrid) .EQ. -1)) THEN |
729 |
IF (indu .GT. 0) THEN |
730 |
C Use the indu value |
731 |
udo(i) = indu - 1 |
732 |
ELSEIF (indu .EQ. -1) THEN |
733 |
C Append one to the current unlimited dim size |
734 |
CALL MNC_DIM_UNLIM_SIZE( fname, unlim_sz, myThid) |
735 |
udo(i) = unlim_sz |
736 |
ELSE |
737 |
C Use the current unlimited dim size |
738 |
CALL MNC_DIM_UNLIM_SIZE( fname, unlim_sz, myThid) |
739 |
udo(i) = unlim_sz - 1 |
740 |
ENDIF |
741 |
ENDIF |
742 |
ENDDO |
743 |
IF (bidim .GT. 0) THEN |
744 |
s(bidim) = lbi |
745 |
e(bidim) = lbi |
746 |
ENDIF |
747 |
IF (bjdim .GT. 0) THEN |
748 |
s(bjdim) = lbj |
749 |
e(bjdim) = lbj |
750 |
ENDIF |
751 |
|
752 |
C DO i = 9,1,-1 |
753 |
C write(*,*) 'i,p(i),s(i),e(i) = ', i,': ',p(i),s(i),e(i) |
754 |
C ENDDO |
755 |
|
756 |
CALL MNC_FILE_ENDDEF(fname, myThid) |
757 |
|
758 |
write(msgbuf,'(5a)') 'reading variable type ''', |
759 |
& vtype(nvf:nvl), ''' within file ''', |
760 |
& fname(1:nfname), '''' |
761 |
|
762 |
C Read the variable one vector at a time |
763 |
DO j7 = s(7),e(7) |
764 |
k7 = (j7 - 1)*p(6) |
765 |
vstart(7) = udo(7) + j7 - s(7) + 1 |
766 |
vcount(7) = 1 |
767 |
DO j6 = s(6),e(6) |
768 |
k6 = (j6 - 1)*p(5) + k7 |
769 |
vstart(6) = udo(6) + j6 - s(6) + 1 |
770 |
vcount(6) = 1 |
771 |
DO j5 = s(5),e(5) |
772 |
k5 = (j5 - 1)*p(4) + k6 |
773 |
vstart(5) = udo(5) + j5 - s(5) + 1 |
774 |
vcount(5) = 1 |
775 |
DO j4 = s(4),e(4) |
776 |
k4 = (j4 - 1)*p(3) + k5 |
777 |
vstart(4) = udo(4) + j4 - s(4) + 1 |
778 |
vcount(4) = 1 |
779 |
DO j3 = s(3),e(3) |
780 |
k3 = (j3 - 1)*p(2) + k4 |
781 |
vstart(3) = udo(3) + j3 - s(3) + 1 |
782 |
vcount(3) = 1 |
783 |
DO j2 = s(2),e(2) |
784 |
k2 = (j2 - 1)*p(1) + k3 |
785 |
vstart(2) = udo(2) + j2 - s(2) + 1 |
786 |
vcount(2) = 1 |
787 |
|
788 |
kr = 0 |
789 |
vstart(1) = udo(1) + 1 |
790 |
vcount(1) = e(1) - s(1) + 1 |
791 |
|
792 |
IF (stype(1:1) .EQ. 'D') THEN |
793 |
err = NF_GET_VARA_DOUBLE(fid, idv, vstart, vcount, resh_d) |
794 |
CALL MNC_HANDLE_ERR(err, msgbuf, myThid) |
795 |
DO j1 = s(1),e(1) |
796 |
k1 = k2 + j1 |
797 |
kr = kr + 1 |
798 |
var(k1) = resh_d(kr) |
799 |
ENDDO |
800 |
ENDIF |
801 |
IF (stype(1:1) .EQ. 'R') THEN |
802 |
err = NF_GET_VARA_REAL(fid, idv, vstart, vcount, resh_r) |
803 |
CALL MNC_HANDLE_ERR(err, msgbuf, myThid) |
804 |
DO j1 = s(1),e(1) |
805 |
k1 = k2 + j1 |
806 |
kr = kr + 1 |
807 |
var(k1) = resh_r(kr) |
808 |
ENDDO |
809 |
ENDIF |
810 |
IF (stype(1:1) .EQ. 'I') THEN |
811 |
err = NF_GET_VARA_INT(fid, idv, vstart, vcount, resh_i) |
812 |
CALL MNC_HANDLE_ERR(err, msgbuf, myThid) |
813 |
DO j1 = s(1),e(1) |
814 |
k1 = k2 + j1 |
815 |
kr = kr + 1 |
816 |
var(k1) = resh_i(kr) |
817 |
ENDDO |
818 |
ENDIF |
819 |
|
820 |
|
821 |
|
822 |
ENDDO |
823 |
ENDDO |
824 |
ENDDO |
825 |
ENDDO |
826 |
ENDDO |
827 |
ENDDO |
828 |
|
829 |
C Close the file |
830 |
CALL MNC_FILE_CLOSE(fname, myThid) |
831 |
|
832 |
C End the lbj,lbi loops |
833 |
ENDDO |
834 |
ENDDO |
835 |
|
836 |
_END_MASTER( myThid ) |
837 |
|
838 |
RETURN |
839 |
END |
840 |
|
841 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
842 |
|
843 |
CEH3 ;;; Local Variables: *** |
844 |
CEH3 ;;; mode:fortran *** |
845 |
CEH3 ;;; End: *** |