/[MITgcm]/MITgcm/pkg/fizhi/fizhi_readwrite_vegtiles.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/fizhi_readwrite_vegtiles.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.20 - (show annotations) (download)
Sun Jun 28 01:05:41 2009 UTC (14 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61s, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.19: +5 -5 lines
add bj in exch2 arrays and S/R

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_readwrite_vegtiles.F,v 1.19 2009/05/12 19:56:35 jmc Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
6 CBOP
7 C !ROUTINE: FIZHI_WRITE_VEGTILES
8 C !INTERFACE:
9 SUBROUTINE FIZHI_WRITE_VEGTILES(fn,pickupflg,myTime,myIter,myThid)
10
11 C !DESCRIPTION:
12
13 C !USES:
14 IMPLICIT NONE
15 #include "SIZE.h"
16 #include "fizhi_SIZE.h"
17 #include "fizhi_land_SIZE.h"
18 #include "fizhi_coms.h"
19 #include "fizhi_land_coms.h"
20 #include "fizhi_earth_coms.h"
21 #include "EEPARAMS.h"
22 #include "PARAMS.h"
23 #ifdef ALLOW_MNC
24 #include "MNC_PARAMS.h"
25 #endif
26 #ifdef ALLOW_EXCH2
27 #include "W2_EXCH2_SIZE.h"
28 #include "W2_EXCH2_TOPOLOGY.h"
29 #endif /* ALLOW_EXCH2 */
30
31 EXTERNAL ILNBLNK
32 INTEGER ILNBLNK
33 INTEGER MDS_RECLEN
34
35 C !INPUT/OUTPUT PARAMETERS:
36 CHARACTER*(*) fn
37 INTEGER pickupflg
38 _RL myTime
39 INTEGER myIter
40 INTEGER myThid
41
42 CEOP
43 C !LOCAL VARIABLES:
44 CHARACTER*1 prec
45 CHARACTER*80 bnam
46 character*(80) dataFName
47 integer ilst
48 integer i,k,n
49 integer ig,jg,tn,iunit
50 integer length_of_rec
51 integer bi,bj,irec,fileprec
52 Real*8 r8seg(nchp)
53
54 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
55
56 DO i = 1,80
57 bnam(i:i) = ' '
58 ENDDO
59 ilst = ILNBLNK(fn)
60 if (pickupflg.eq.0) then
61 prec = 'D'
62 fileprec = 64
63 WRITE(bnam,'(a,a)') 'pickup_vegtiles.', fn(1:ilst)
64 else
65 prec = 'D'
66 fileprec = 64
67 WRITE(bnam,'(a,a)') 'state_vegtiles.', fn(1:ilst)
68 endif
69
70 #ifdef ALLOW_MNC
71 IF (useMNC.AND. pickup_write_mnc) THEN
72
73 C Write fizhi veg-space variables using the MNC package
74 CALL MNC_CW_SET_UDIM(bnam, 1, myThid)
75 CALL MNC_CW_RL_W_S('D',bnam,0,0,'T', myTime, myThid)
76 CALL MNC_CW_I_W_S('I',bnam,0,0,'iter',myIter,myThid)
77
78 C fizhi_coms.h
79 CALL MNC_CW_RL_W(prec,bnam,0,0,'ctmt', ctmt, myThid)
80 CALL MNC_CW_RL_W(prec,bnam,0,0,'xxmt', xxmt, myThid)
81 CALL MNC_CW_RL_W(prec,bnam,0,0,'yymt', yymt, myThid)
82 CALL MNC_CW_RL_W(prec,bnam,0,0,'zetamt', zetamt, myThid)
83 CALL MNC_CW_RL_W(prec,bnam,0,0,'xlmt', xlmt, myThid)
84 CALL MNC_CW_RL_W(prec,bnam,0,0,'khmt', khmt, myThid)
85 CALL MNC_CW_RL_W(prec,bnam,0,0,'tke', tke, myThid)
86
87 C fizhi_land_coms.h
88 CALL MNC_CW_RL_W(prec,bnam,0,0,'tcanopy', tcanopy, myThid)
89 CALL MNC_CW_RL_W(prec,bnam,0,0,'tdeep', tdeep, myThid)
90 CALL MNC_CW_RL_W(prec,bnam,0,0,'ecanopy', ecanopy, myThid)
91 CALL MNC_CW_RL_W(prec,bnam,0,0,'swetshal', swetshal, myThid)
92 CALL MNC_CW_RL_W(prec,bnam,0,0,'swetroot', swetroot, myThid)
93 CALL MNC_CW_RL_W(prec,bnam,0,0,'swetdeep', swetdeep, myThid)
94 CALL MNC_CW_RL_W(prec,bnam,0,0,'snodep', snodep, myThid)
95 CALL MNC_CW_RL_W(prec,bnam,0,0,'capac', capac, myThid)
96 CALL MNC_CW_RL_W(prec,bnam,0,0,'chlt', chlt, myThid)
97 CALL MNC_CW_RL_W(prec,bnam,0,0,'chlon', chlon, myThid)
98 CALL MNC_CW_I_W('I',bnam,0,0,'igrd', igrd, myThid)
99
100 C fizhi_earth_coms.h
101 CALL MNC_CW_I_W('I',bnam,0,0,'ityp', ityp, myThid)
102 CALL MNC_CW_RL_W(prec,bnam,0,0,'chfr', chfr, myThid)
103
104 ENDIF
105 #endif /* Not ALLOW_MNC sequence */
106
107
108 call MDSFINDUNIT( iunit, mythid )
109 length_of_rec=MDS_RECLEN( fileprec, nchp, mythid )
110
111 DO bj = myByLo(myThid), myByHi(myThid)
112 DO bi = myBxLo(myThid), myBxHi(myThid)
113
114 #ifdef ALLOW_EXCH2
115 tn = W2_myTileList(bi,bj)
116 iG = tn
117 jG = 1
118 #else
119 iG = bi+(myXGlobalLo-1)/sNx
120 jG = bj+(myYGlobalLo-1)/sNy
121 tn = (jG - 1)*(nPx*nSx) + iG
122 #endif /* ALLOW_EXCH2 */
123
124 write(dataFname(1:80),'(a,2a,i3.3,a,i3.3,a)')
125 & 'pickup_vegtiles.',fn(1:ilst),'.',iG,'.',jG,'.data'
126 open( iUnit, file=dataFName, status='unknown',
127 & access='direct', recl=length_of_rec )
128
129 C First write single-level turbulence fields
130 do n = 1,nchp
131 r8seg(n) = ctmt(n,bi,bj)
132 enddo
133 #ifdef _BYTESWAPIO
134 call MDS_BYTESWAPR8( nchp, r8seg )
135 #endif
136 write(iunit,rec=1) r8seg
137
138 do n = 1,nchp
139 r8seg(n) = xxmt(n,bi,bj)
140 enddo
141 #ifdef _BYTESWAPIO
142 call MDS_BYTESWAPR8( nchp, r8seg )
143 #endif
144 write(iunit,rec=2) r8seg
145
146 do n = 1,nchp
147 r8seg(n) = yymt(n,bi,bj)
148 enddo
149 #ifdef _BYTESWAPIO
150 call MDS_BYTESWAPR8( nchp, r8seg )
151 #endif
152 write(iunit,rec=3) r8seg
153
154 do n = 1,nchp
155 r8seg(n) = zetamt(n,bi,bj)
156 enddo
157 #ifdef _BYTESWAPIO
158 call MDS_BYTESWAPR8( nchp, r8seg )
159 #endif
160 write(iunit,rec=4) r8seg
161
162 C And now write Multi-level turbulence fields
163 do k = 1,Nrphys
164 do n = 1,nchp
165 r8seg(n) = xlmt(n,k,bi,bj)
166 enddo
167 #ifdef _BYTESWAPIO
168 call MDS_BYTESWAPR8( nchp, r8seg )
169 #endif
170 irec = 4 + 0*Nrphys + k
171 write(iunit,rec=irec) r8seg
172 enddo
173
174 do k = 1,Nrphys
175 do n = 1,nchp
176 r8seg(n) = khmt(n,k,bi,bj)
177 enddo
178 #ifdef _BYTESWAPIO
179 call MDS_BYTESWAPR8( nchp, r8seg )
180 #endif
181 irec = 4 + 1*Nrphys + k
182 write(iunit,rec=irec) r8seg
183 enddo
184
185 do k = 1,Nrphys
186 do n = 1,nchp
187 r8seg(n) = tke(n,k,bi,bj)
188 enddo
189 #ifdef _BYTESWAPIO
190 call MDS_BYTESWAPR8( nchp, r8seg )
191 #endif
192 irec = 4 + 2*Nrphys + k
193 write(iunit,rec=irec) r8seg
194 enddo
195
196 C And finally, write land surface fields
197 do n = 1,nchp
198 r8seg(n) = tcanopy(n,bi,bj)
199 enddo
200 #ifdef _BYTESWAPIO
201 call MDS_BYTESWAPR8( nchp, r8seg )
202 #endif
203 irec = 4 + 3*Nrphys + 1
204 write(iunit,rec=irec) r8seg
205
206 do n = 1,nchp
207 r8seg(n) = tdeep(n,bi,bj)
208 enddo
209 #ifdef _BYTESWAPIO
210 call MDS_BYTESWAPR8( nchp, r8seg )
211 #endif
212 irec = 4 + 3*Nrphys + 2
213 write(iunit,rec=irec) r8seg
214
215 do n = 1,nchp
216 r8seg(n) = ecanopy(n,bi,bj)
217 enddo
218 #ifdef _BYTESWAPIO
219 call MDS_BYTESWAPR8( nchp, r8seg )
220 #endif
221 irec = 4 + 3*Nrphys + 3
222 write(iunit,rec=irec) r8seg
223
224 do n = 1,nchp
225 r8seg(n) = swetshal(n,bi,bj)
226 enddo
227 #ifdef _BYTESWAPIO
228 call MDS_BYTESWAPR8( nchp, r8seg )
229 #endif
230 irec = 4 + 3*Nrphys + 4
231 write(iunit,rec=irec) r8seg
232
233 do n = 1,nchp
234 r8seg(n) = swetroot(n,bi,bj)
235 enddo
236 #ifdef _BYTESWAPIO
237 call MDS_BYTESWAPR8( nchp, r8seg )
238 #endif
239 irec = 4 + 3*Nrphys + 5
240 write(iunit,rec=irec) r8seg
241
242 do n = 1,nchp
243 r8seg(n) = swetdeep(n,bi,bj)
244 enddo
245 #ifdef _BYTESWAPIO
246 call MDS_BYTESWAPR8( nchp, r8seg )
247 #endif
248 irec = 4 + 3*Nrphys + 6
249 write(iunit,rec=irec) r8seg
250
251 do n = 1,nchp
252 r8seg(n) = snodep(n,bi,bj)
253 enddo
254 #ifdef _BYTESWAPIO
255 call MDS_BYTESWAPR8( nchp, r8seg )
256 #endif
257 irec = 4 + 3*Nrphys + 7
258 write(iunit,rec=irec) r8seg
259
260 do n = 1,nchp
261 r8seg(n) = capac(n,bi,bj)
262 enddo
263 #ifdef _BYTESWAPIO
264 call MDS_BYTESWAPR8( nchp, r8seg )
265 #endif
266 irec = 4 + 3*Nrphys + 8
267 write(iunit,rec=irec) r8seg
268
269 close(iunit)
270
271 C End of bi bj loop
272 enddo
273 enddo
274
275 RETURN
276 END
277
278
279 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
280 CBOP
281 C !ROUTINE: FIZHI_READ_VEGTILES
282 C !INTERFACE:
283 SUBROUTINE FIZHI_READ_VEGTILES(Iter,prec,myThid)
284
285 C !DESCRIPTION:
286
287 C !USES:
288 IMPLICIT NONE
289 #include "SIZE.h"
290 #include "fizhi_SIZE.h"
291 #include "fizhi_land_SIZE.h"
292 #include "fizhi_coms.h"
293 #include "fizhi_land_coms.h"
294 #include "fizhi_earth_coms.h"
295 #include "EEPARAMS.h"
296 #include "PARAMS.h"
297 #ifdef ALLOW_MNC
298 #include "MNC_PARAMS.h"
299 #endif
300 #ifdef ALLOW_EXCH2
301 #include "W2_EXCH2_SIZE.h"
302 #include "W2_EXCH2_TOPOLOGY.h"
303 #endif /* ALLOW_EXCH2 */
304
305 EXTERNAL ILNBLNK
306 INTEGER ILNBLNK
307 INTEGER MDS_RECLEN
308
309 C !INPUT/OUTPUT PARAMETERS:
310 CHARACTER*1 prec
311 INTEGER Iter
312 INTEGER myThid
313
314 CEOP
315 C !LOCAL VARIABLES:
316 CHARACTER*80 fn
317 CHARACTER*80 bnam
318 integer ilst
319 character*(80) dataFName
320 integer i,k,n
321 integer ig,jg,tn,iunit
322 integer length_of_rec
323 integer bi,bj,irec,fileprec
324 Real*8 r8seg(nchp)
325
326 DO i = 1,80
327 bnam(i:i) = ' '
328 ENDDO
329 WRITE(fn,'(a,I10.10)') 'pickup_vegtiles.',Iter
330 ilst = ILNBLNK(fn)
331 WRITE(bnam,'(a,I10.10)') 'pickup_vegtiles.',Iter
332 fileprec = 64
333
334 #ifdef ALLOW_MNC
335 IF (useMNC.AND. pickup_write_mnc) THEN
336
337 C Write fizhi veg-space variables using the MNC package
338 CALL MNC_FILE_CLOSE_ALL_MATCHING(bnam, myThid)
339 CALL MNC_CW_SET_UDIM(bnam, 1, myThid)
340
341 C fizhi_coms.h
342 CALL MNC_CW_RL_R(prec,bnam,0,0,'ctmt', ctmt, myThid)
343 CALL MNC_CW_RL_R(prec,bnam,0,0,'xxmt', xxmt, myThid)
344 CALL MNC_CW_RL_R(prec,bnam,0,0,'yymt', yymt, myThid)
345 CALL MNC_CW_RL_R(prec,bnam,0,0,'zetamt', zetamt, myThid)
346 CALL MNC_CW_RL_R(prec,bnam,0,0,'xlmt', xlmt, myThid)
347 CALL MNC_CW_RL_R(prec,bnam,0,0,'khmt', khmt, myThid)
348 CALL MNC_CW_RL_R(prec,bnam,0,0,'tke', tke, myThid)
349
350 C fizhi_land_coms.h
351 CALL MNC_CW_RL_R(prec,bnam,0,0,'tcanopy', tcanopy, myThid)
352 CALL MNC_CW_RL_R(prec,bnam,0,0,'tdeep', tdeep, myThid)
353 CALL MNC_CW_RL_R(prec,bnam,0,0,'ecanopy', ecanopy, myThid)
354 CALL MNC_CW_RL_R(prec,bnam,0,0,'swetshal', swetshal, myThid)
355 CALL MNC_CW_RL_R(prec,bnam,0,0,'swetroot', swetroot, myThid)
356 CALL MNC_CW_RL_R(prec,bnam,0,0,'swetdeep', swetdeep, myThid)
357 CALL MNC_CW_RL_R(prec,bnam,0,0,'snodep', snodep, myThid)
358 CALL MNC_CW_RL_R(prec,bnam,0,0,'capac', capac, myThid)
359
360 ENDIF
361 #endif /* Not ALLOW_MNC sequence */
362
363 call MDSFINDUNIT( iunit, mythid )
364 length_of_rec=MDS_RECLEN( fileprec, nchp, mythid )
365
366 DO bj = myByLo(myThid), myByHi(myThid)
367 DO bi = myBxLo(myThid), myBxHi(myThid)
368
369 #ifdef ALLOW_EXCH2
370 tn = W2_myTileList(bi,bj)
371 iG = tn
372 jG = 1
373 #else
374 iG = bi+(myXGlobalLo-1)/sNx
375 jG = bj+(myYGlobalLo-1)/sNy
376 tn = (jG - 1)*(nPx*nSx) + iG
377 #endif /* ALLOW_EXCH2 */
378
379 write(dataFname(1:80),'(2a,i3.3,a,i3.3,a)')
380 & fn(1:ilst),'.',iG,'.',jG,'.data'
381 print *,' Opening ',dataFName
382 open( iUnit, file=dataFName, status='old',
383 & access='direct', recl=length_of_rec )
384
385 irec = 0
386 C First read single-level turbulence fields
387 read(iunit,rec=1) r8seg
388 #ifdef _BYTESWAPIO
389 call MDS_BYTESWAPR8( nchp, r8seg )
390 #endif
391 do n = 1,nchp
392 ctmt(n,bi,bj) = r8seg(n)
393 enddo
394
395 read(iunit,rec=2) r8seg
396 #ifdef _BYTESWAPIO
397 call MDS_BYTESWAPR8( nchp, r8seg )
398 #endif
399 do n = 1,nchp
400 xxmt(n,bi,bj) = r8seg(n)
401 enddo
402
403 read(iunit,rec=3) r8seg
404 #ifdef _BYTESWAPIO
405 call MDS_BYTESWAPR8( nchp, r8seg )
406 #endif
407 do n = 1,nchp
408 yymt(n,bi,bj) = r8seg(n)
409 enddo
410
411 read(iunit,rec=4) r8seg
412 #ifdef _BYTESWAPIO
413 call MDS_BYTESWAPR8( nchp, r8seg )
414 #endif
415 do n = 1,nchp
416 zetamt(n,bi,bj) = r8seg(n)
417 enddo
418
419 C And now read Multi-level turbulence fields
420 do k = 1,Nrphys
421 irec = 4 + 0*Nrphys + k
422 read(iunit,rec=irec) r8seg
423 #ifdef _BYTESWAPIO
424 call MDS_BYTESWAPR8( nchp, r8seg )
425 #endif
426 do n = 1,nchp
427 xlmt(n,k,bi,bj) = r8seg(n)
428 enddo
429 enddo
430
431 do k = 1,Nrphys
432 irec = 4 + 1*Nrphys + k
433 read(iunit,rec=irec) r8seg
434 #ifdef _BYTESWAPIO
435 call MDS_BYTESWAPR8( nchp, r8seg )
436 #endif
437 do n = 1,nchp
438 khmt(n,k,bi,bj) = r8seg(n)
439 enddo
440 enddo
441
442 do k = 1,Nrphys
443 irec = 4 + 2*Nrphys + k
444 read(iunit,rec=irec) r8seg
445 #ifdef _BYTESWAPIO
446 call MDS_BYTESWAPR8( nchp, r8seg )
447 #endif
448 do n = 1,nchp
449 tke(n,k,bi,bj) = r8seg(n)
450 enddo
451 enddo
452
453 C And finally, read land surface fields
454 irec = 4 + 3*Nrphys + 1
455 read(iunit,rec=irec) r8seg
456 #ifdef _BYTESWAPIO
457 call MDS_BYTESWAPR8( nchp, r8seg )
458 #endif
459 do n = 1,nchp
460 tcanopy(n,bi,bj) = r8seg(n)
461 enddo
462
463 irec = 4 + 3*Nrphys + 2
464 read(iunit,rec=irec) r8seg
465 #ifdef _BYTESWAPIO
466 call MDS_BYTESWAPR8( nchp, r8seg )
467 #endif
468 do n = 1,nchp
469 tdeep(n,bi,bj) = r8seg(n)
470 enddo
471
472 irec = 4 + 3*Nrphys + 3
473 read(iunit,rec=irec) r8seg
474 #ifdef _BYTESWAPIO
475 call MDS_BYTESWAPR8( nchp, r8seg )
476 #endif
477 do n = 1,nchp
478 ecanopy(n,bi,bj) = r8seg(n)
479 enddo
480
481 irec = 4 + 3*Nrphys + 4
482 read(iunit,rec=irec) r8seg
483 #ifdef _BYTESWAPIO
484 call MDS_BYTESWAPR8( nchp, r8seg )
485 #endif
486 do n = 1,nchp
487 swetshal(n,bi,bj) = r8seg(n)
488 enddo
489
490 irec = 4 + 3*Nrphys + 5
491 read(iunit,rec=irec) r8seg
492 #ifdef _BYTESWAPIO
493 call MDS_BYTESWAPR8( nchp, r8seg )
494 #endif
495 do n = 1,nchp
496 swetroot(n,bi,bj) = r8seg(n)
497 enddo
498
499 irec = 4 + 3*Nrphys + 6
500 read(iunit,rec=irec) r8seg
501 #ifdef _BYTESWAPIO
502 call MDS_BYTESWAPR8( nchp, r8seg )
503 #endif
504 do n = 1,nchp
505 swetdeep(n,bi,bj) = r8seg(n)
506 enddo
507
508 irec = 4 + 3*Nrphys + 7
509 read(iunit,rec=irec) r8seg
510 #ifdef _BYTESWAPIO
511 call MDS_BYTESWAPR8( nchp, r8seg )
512 #endif
513 do n = 1,nchp
514 snodep(n,bi,bj) = r8seg(n)
515 enddo
516
517 irec = 4 + 3*Nrphys + 8
518 read(iunit,rec=irec) r8seg
519 #ifdef _BYTESWAPIO
520 call MDS_BYTESWAPR8( nchp, r8seg )
521 #endif
522 do n = 1,nchp
523 capac(n,bi,bj) = r8seg(n)
524 enddo
525
526 close(iunit)
527
528 C End of bi bj loop
529 enddo
530 enddo
531
532
533 RETURN
534 END
535
536 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|

  ViewVC Help
Powered by ViewVC 1.1.22