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