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