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