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

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

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


Revision 1.20 - (show annotations) (download)
Tue Mar 8 22:14:20 2005 UTC (19 years, 3 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint57g_post, checkpoint57g_pre, checkpoint57f_pre, checkpoint57f_post
Changes since 1.19: +104 -88 lines
Read new sst files

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/update_ocean_exports.F,v 1.19 2005/03/08 20:05:46 molod Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 subroutine update_ocean_exports (myTime, myIter, myThid)
6 c----------------------------------------------------------------------
7 c Subroutine update_ocean_exports - 'Wrapper' routine to update
8 c the fields related to the ocean's surface that are needed
9 c by fizhi (sst and sea ice extent).
10 c
11 c Call: getsst (Return the current sst field-read dataset if needed)
12 c getsice (Return the current sea ice field-read data if needed)
13 c-----------------------------------------------------------------------
14 implicit none
15 #include "SIZE.h"
16 #include "GRID.h"
17 #include "fizhi_ocean_coms.h"
18 #include "EEPARAMS.h"
19 #include "chronos.h"
20 #ifdef ALLOW_EXCH2
21 #include "W2_EXCH2_TOPOLOGY.h"
22 #include "W2_EXCH2_PARAMS.h"
23 #endif /* ALLOW_EXCH2 */
24
25 integer myIter, myThid
26 _RL myTime
27
28 integer i, j, bi, bj, bislot, bjslot
29 integer im1, im2, jm1, jm2, idim1, idim2, jdim1, jdim2
30 integer xsize, ysize
31 integer ksst,kice
32 _RL sstmin
33 parameter ( sstmin = 273.16 )
34
35 #if defined(ALLOW_EXCH2)
36 PARAMETER ( xsize = exch2_domain_nxt * sNx )
37 PARAMETER ( ysize = exch2_domain_nyt * sNy )
38 #else
39 PARAMETER ( xsize = Nx )
40 PARAMETER ( ysize = Ny )
41 #endif
42 _RL sst1(xsize,ysize),sst2(xsize,ysize)
43 _RL sice1(xsize,ysize),sice2(xsize,ysize)
44
45 save sst1, sst2, sice1, sice2
46
47 idim1 = 1-OLx
48 idim2 = sNx+OLx
49 jdim1 = 1-OLy
50 jdim2 = sNy+OLy
51 im1 = 1
52 im2 = sNx
53 jm1 = 1
54 jm2 = sNy
55
56 call mdsfindunit( ksst, myThid )
57 call mdsfindunit( kice, myThid )
58
59 C***********************************************************************
60
61 DO BJ = myByLo(myThid),myByHi(myThid)
62 DO BI = myBxLo(myThid),myBxHi(myThid)
63 #if defined(ALLOW_EXCH2)
64 bislot = exch2_txglobalo(W2_myTileList(bi))-1
65 bjslot = exch2_tyglobalo(W2_myTileList(bi))-1
66 #else
67 bislot = myXGlobalLo-1+(bi-1)*sNx
68 bjslot = myYGlobalLo-1+(bj-1)*sNy
69 #endif
70
71 call getsst(ksst,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,nSx,
72 . nSy,xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms,sst1,sst2,sst)
73 call getsice(kice,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,nSx,
74 . nSy,xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms,sice1,sice2,sice)
75
76 c Check for Minimum Open-Water SST
77 c --------------------------------
78 do j=jm1,jm2
79 do i=im1,im2
80 if(sice(i,j,bi,bj).eq.0.0 .and. sst(i,j,bi,bj).lt.sstmin)
81 . sst(i,j,bi,bj) = sstmin
82 enddo
83 enddo
84
85 ENDDO
86 ENDDO
87 _EXCH_XYZ_R8(sst,myThid)
88 _EXCH_XYZ_R8(sice,myThid)
89
90 return
91 end
92
93 subroutine getsice(iunit,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,
94 . nSumx,nSumy,xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms,
95 . sicebc1,sicebc2,sice)
96 C***********************************************************************
97 C
98 C!ROUTINE: GETSICE
99 C!DESCRIPTION: GETSICE returns the sea ice depth.
100 C! This routine is adaptable for any frequency
101 C! data upto a daily frequency.
102 C! note: for diurnal data ndmax should be increased.
103 C
104 C!INPUT PARAMETERS:
105 C! iunit Unit number assigned to the sice data file
106 C! idim1 Start dimension in x-direction
107 C! idim2 End dimension in x-direction
108 C! jdim1 Start dimension in y-direction
109 C! jdim2 End dimension in y-direction
110 C! im1 Begin of x-direction span for filling sice
111 C! im2 End of x-direction span for filling sice
112 C! jm1 Begin of y-direction span for filling sice
113 C! jm2 End of y-direction span for filling sice
114 C! nSumx Number of processors in x-direction (local processor)
115 C! nSumy Number of processors in y-direction (local processor)
116 C! xsize Number of processors in x-direction (global)
117 C! ysize Number of processors in y-direction (global)
118 C! bi Processor number in x-direction (local to processor)
119 C! bj Processor number in y-direction (local to processor)
120 C! bislot Processor number in x-direction (global)
121 C! bjslot Processor number in y-direction (global)
122 C! nymd YYMMDD of the current model timestep
123 C! nhms HHMMSS of the model time
124 C
125 C!OUTPUT PARAMETERS:
126 C! sice(idim1:idim2,jdim1:jdim2,nSumx,nSumy) Sea ice depth in meters
127 C
128 C!ROUTINES CALLED:
129 C
130 C! bcdata Reads the data for a given unit number
131 C! bcheader Reads the header info for a given unit number
132 C! interp_time Returns weights for linear interpolation
133 C
134 C--------------------------------------------------------------------------
135
136 implicit none
137 #include "SIZE.h"
138
139 integer iunit,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,nSumx,nSumy
140 integer xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms
141
142 _RL sicebc1(xsize,ysize)
143 _RL sicebc2(xsize,ysize)
144 _RL sice(idim1:idim2,jdim1:jdim2,nSumx,nSumy)
145
146 C Maximum number of dates in one year for the data
147 integer ndmax
148 parameter (ndmax = 370)
149
150 character*8 cname
151 character*80 cdscrip
152 character*40 sicedata
153 _RL fac1, fac2, lat0, lon0, timebc1, timebc2, timemod, undef
154 logical first, found, error
155 integer i,j,n,nn,iyear,iyearbc,nd,ndby3,imbc,jmbc
156 integer ndatebc,nhmsbc(ndmax), nhmsbc1, nhmsbc2,nrec
157 integer nymdbc(ndmax),nymdbc1,nymdbc2,nymdmod
158
159
160 C--------- Variable Initialization ---------------------------------
161
162 data first /.true./
163 data error /.false./
164
165 c save header info
166 save imbc,jmbc,lat0,lon0,ndatebc,undef,nymdbc,nhmsbc
167 save nymdbc1, nymdbc2, nhmsbc1, nhmsbc2
168 save first
169
170 c this only works for between 1950-2050
171 if (nymd .lt. 500101) then
172 nymdmod = 20000000 + nymd
173 else if (nymd .le. 991231) then
174 nymdmod = 19000000 + nymd
175 else
176 nymdmod = nymd
177 endif
178
179 sicedata='sice19232.weekly.clim'
180
181 c initialize so that first time through they have values for the check
182 c these vaules make the iyear .ne. iyearbc true anyways for
183 c for the first time so first isnt checked below.
184
185 if (first) then
186 nymdbc(2) = 0
187 nymdbc1 = 0
188 nymdbc2 = 0
189 nhmsbc1 = 0
190 nhmsbc2 = 0
191 c first = .false.
192 endif
193
194 C---------- Read in Header file ----------------------------------
195
196 iyear = nymdmod/10000
197 iyearbc = nymdbc(2)/10000
198
199 if( iyear.ne.iyearbc ) then
200
201 close(iunit)
202 open (iunit,file=sicedata,form='unformatted',access='direct',
203 . recl=xsize*ysize*4)
204 nrec = 1
205 call bcheader (iunit, ndmax, nrec,
206 . cname, cdscrip, imbc, jmbc, lat0, lon0,
207 . ndatebc, nymdbc, nhmsbc, undef, error)
208
209 C--------- Check data for Compatibility ------------------------------
210
211 C Check for correct data in boundary condition file
212 if (.not.error .and. cname.ne.'SICE') then
213 write(6,*)'Wrong data in SICE boundary condition file => ',cname
214 error = .true.
215 endif
216
217 C Check Horizontal Resolution
218 if(.not.error.and.imbc*jmbc.ne.xsize*ysize)then
219 write(6,*) ' B.C. Resolution DOES NOT match Model Resolution!'
220 write(6,*) ' B.C. Resolution: ',imbc*jmbc
221 write(6,*) 'Model Resolution: ',xsize*ysize
222 error = .true.
223 endif
224
225 C Check Year
226 iyearbc = nymdbc(2)/10000
227 if (.not.error .and. iyear.ne.iyearbc .and. iyearbc.ne.0) then
228 write(6,*)' B.C. Year DOES NOT match REQUESTED Year!'
229 write(6,*)' B.C. Year: ', iyearbc
230 write(6,*)'Requested Year: ', iyear
231 error = .true.
232 endif
233
234 if (.not.error) then
235 C if climatology, fill dates for data with current model year
236 if (iyearbc.eq.0) then
237 write(6,*)
238 write(6,*) 'Climatological Dataset is being used.'
239 write(6,*) 'Current model year to be used to fill Header Dates'
240 do n = 2, ndatebc-1
241 nymdbc(n) = nymdbc(n) +(nymdmod/10000)*10000
242 enddo
243 C For the first date subtract 1 year from the current model NYMD
244 n = 1
245 nymdbc(n) = nymdbc(n) +(nymdmod/10000-1)*10000
246 C For the last date add 1 year to the current model NYMD
247 n = ndatebc
248 nymdbc(n) = nymdbc(n) +(nymdmod/10000+1)*10000
249 endif
250
251 C Write out header info
252 write(6,*) ' Updated boundary condition data'
253 write(6,*) ' ---------------------------------'
254 write(6,*) ' Variable: ',cname
255 write(6,*) ' Description: ',cdscrip
256 write(6,*) ' Resolution: x= ',imbc,' y= ',jmbc,
257 . ' Undefined value = ',undef
258 write(6,*) ' Starting latitude = ',lat0,' longitude =',lon0
259 write(6,*) ' Data valid at these times: '
260 ndby3 = ndatebc/3
261 do n = 1, ndby3*3,3
262 write(6,1000) (nn,nymdbc(nn),nhmsbc(nn),nn=n,n+2)
263 1000 format(3(2x,i3,':',i8,2x,i8))
264 enddo
265 write(6,1000) (nn,nymdbc(nn),nhmsbc(nn),nn=ndby3*3+1,ndatebc)
266 endif
267
268 endif
269
270 C---------- Read sice data if necessary -------------------------------
271
272 found = .false.
273 nd = 2
274
275 c If model time is not within the times of saved sice data
276 c from previous call to getsice then read new data
277
278 timemod = float(nymdmod) + float(nhms) /1000000
279 timebc1 = float(nymdbc1) + float(nhmsbc1)/1000000
280 timebc2 = float(nymdbc2) + float(nhmsbc2)/1000000
281
282 if (timemod .lt. timebc1 .or. timemod .ge. timebc2) then
283
284 do while (.not.found .and. nd .le. ndatebc)
285 timebc2 = float(nymdbc(nd)) + float(nhmsbc(nd))/1000000
286 if (timebc2 .gt. timemod) then
287 nymdbc1 = nymdbc(nd-1)
288 nymdbc2 = nymdbc(nd)
289 nhmsbc1 = nhmsbc(nd-1)
290 nhmsbc2 = nhmsbc(nd)
291 call bcdata (iunit,imbc,jmbc,nd,nd+1,sicebc1,sicebc2)
292 found = .true.
293 else
294 nd = nd + 1
295 endif
296 enddo
297
298 c Otherwise the data from the last time in getsice surrounds the
299 c current model time.
300
301 else
302 found = .true.
303 endif
304
305 if (.not.found) then
306 print *, 'STOP: Could not find SICE dates for model time.'
307 call my_finalize
308 call my_exit (101)
309 endif
310
311 C---------- Interpolate sice data ------------------------------------
312
313 call interp_time(nymdmod,nhms,nymdbc1,nhmsbc1,nymdbc2,nhmsbc2,
314 . fac1,fac2)
315
316 do j = jm1,jm2
317 do i = im1,im2
318 sice(i,j,bi,bj) = sicebc1(i+bislot,j+bjslot)*fac1
319 . + sicebc2(i+bislot,j+bjslot)*fac2
320 c average to 0 or 1
321 c -----------------
322 if (sice(i,j,bi,bj) .ge. 0.5) then
323 sice(i,j,bi,bj) = 1.
324 else
325 sice(i,j,bi,bj) = 0.
326 endif
327 enddo
328 enddo
329
330 C---------- Fill sice with depth of ice ------------------------------------
331 do j = jm1,jm2
332 do i = im1,im2
333 if (sice(i,j,bi,bj) .eq. 1.) then
334 sice(i,j,bi,bj) = 3.
335 endif
336 enddo
337 enddo
338 C---------------------------------------------------------------------------
339
340 return
341 end
342 subroutine getsst(iunit,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,
343 . nSumx,nSumy,xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms,
344 . sstbc1,sstbc2,sst)
345 C***********************************************************************
346 C
347 C!ROUTINE: GETSST
348 C!DESCRIPTION: GETSST gets the SST data.
349 C! This routine is adaptable for any frequency
350 C! data upto a daily frequency.
351 C! note: for diurnal data ndmax should be increased.
352 C
353 C!INPUT PARAMETERS:
354 C! iunit Unit number assigned to the sice data file
355 C! idim1 Start dimension in x-direction
356 C! idim2 End dimension in x-direction
357 C! jdim1 Start dimension in y-direction
358 C! jdim2 End dimension in y-direction
359 C! im1 Begin of x-direction span for filling sst
360 C! im2 End of x-direction span for filling sst
361 C! jm1 Begin of y-direction span for filling sst
362 C! jm2 End of y-direction span for filling sst
363 C! nSumx Number of processors in x-direction (local processor)
364 C! nSumy Number of processors in y-direction (local processor)
365 C! xsize x-dimension of global array
366 C! ysize y-dimension of global array
367 C! bi Processor number in x-direction (local to processor)
368 C! bj Processor number in y-direction (local to processor)
369 C! bislot Slot number into global array in x-direction (global)
370 C! bjslot Slot number into global array in y-direction (global)
371 C! nymd YYMMDD of the current model timestep
372 C! nhms HHMMSS of the model time
373 C
374 C!OUTPUT PARAMETERS:
375 C! sst(idim1:idim2,jdim1:jdim2,nSumx,nSumy) Sea surface temperature (K)
376 C
377 C!ROUTINES CALLED:
378 C
379 C! bcdata Reads the data for a given unit number
380 C! bcheader Reads the header info for a given unit number
381 C! interp_time Returns weights for linear interpolation
382 C
383 C--------------------------------------------------------------------------
384
385 implicit none
386 #include "SIZE.h"
387
388 integer iunit,idim1,idim2,jdim1,jdim2,im1,im2,jm1,jm2,nSumx,nSumy
389 integer xsize,ysize,bi,bj,bislot,bjslot,nymd,nhms
390
391 _RL sstbc1(xsize,ysize)
392 _RL sstbc2(xsize,ysize)
393 _RL sst(idim1:idim2,jdim1:jdim2,nSumx,nSumy)
394
395 C Maximum number of dates in one year for the data
396 integer ndmax
397 parameter (ndmax = 370)
398
399 character*8 cname
400 character*80 cdscrip
401 character*20 sstdata
402 _RL fac1, fac2, lat0, lon0, timebc1, timebc2, timemod, undef
403 logical first, found, error
404 integer i,j,n,nn,iyear,iyearbc,nd,ndby3,imbc,jmbc
405 integer ndatebc,nhmsbc(ndmax), nhmsbc1, nhmsbc2,nrec
406 integer nymdbc(ndmax),nymdbc1,nymdbc2,nymdmod
407
408
409 C--------- Variable Initialization ---------------------------------
410
411 data first /.true./
412 data error /.false./
413
414 c save header info
415 save imbc,jmbc,lat0,lon0,ndatebc,undef,nymdbc,nhmsbc
416 save nymdbc1, nymdbc2, nhmsbc1, nhmsbc2
417 save first
418
419 c this only works for between 1950-2050
420 if (nymd .lt. 500101) then
421 nymdmod = 20000000 + nymd
422 else if (nymd .le. 991231) then
423 nymdmod = 19000000 + nymd
424 else
425 nymdmod = nymd
426 endif
427
428 sstdata='sst19232.weekly.clim'
429
430 c initialize so that first time through they have values for the check
431 c these vaules make the iyear .ne. iyearbc true anyways for
432 c for the first time so first isnt checked below.
433 if (first) then
434 nymdbc(2) = 0
435 nymdbc1 = 0
436 nymdbc2 = 0
437 nhmsbc1 = 0
438 nhmsbc2 = 0
439 c first = .false.
440 endif
441
442 C---------- Read in Header file ----------------------------------
443
444 iyear = nymdmod/10000
445 iyearbc = nymdbc(2)/10000
446
447 if( iyear.ne.iyearbc ) then
448
449 close(iunit)
450 open (iunit,file=sstdata,form='unformatted',access='direct',
451 . recl=xsize*ysize*4)
452 nrec = 1
453 call bcheader (iunit, ndmax, nrec,
454 . cname, cdscrip, imbc, jmbc, lat0, lon0,
455 . ndatebc, nymdbc, nhmsbc, undef, error)
456
457 C--------- Check data for Compatibility
458
459 C Check for correct data in boundary condition file
460 if (.not.error .and. cname.ne.'SST') then
461 write(6,*)'Wrong data in SST boundary condition file => ',cname
462 error = .true.
463 endif
464
465 C Check Horizontal Resolution
466 if(.not.error.and.imbc*jmbc.ne.xsize*ysize)then
467 write(6,*) ' B.C. Resolution DOES NOT match Model Resolution!'
468 write(6,*) ' B.C. Resolution: ',imbc*jmbc
469 write(6,*) 'Model Resolution: ',xsize*ysize
470 error = .true.
471 endif
472
473 C Check Year
474 iyearbc = nymdbc(2)/10000
475 if (.not.error .and. iyear.ne.iyearbc .and. iyearbc.ne.0) then
476 write(6,*)' B.C. Year DOES NOT match REQUESTED Year!'
477 write(6,*)' B.C. Year: ', iyearbc
478 write(6,*)'Requested Year: ', iyear
479 error = .true.
480 endif
481
482 if (.not.error) then
483 C if climatology, fill dates for data with current model year
484 if (iyearbc.eq.0) then
485 write(6,*)
486 write(6,*)'Climatological Dataset is being used.'
487 write(6,*)'Current model year is used to fill Header Dates'
488 do n = 2, ndatebc-1
489 nymdbc(n) = nymdbc(n) +(nymdmod/10000)*10000
490 enddo
491 C For the first date subtract 1 year from the current model NYMD
492 n = 1
493 nymdbc(n) = nymdbc(n) +(nymdmod/10000-1)*10000
494 C For the last date add 1 year to the current model NYMD
495 n = ndatebc
496 nymdbc(n) = nymdbc(n) +(nymdmod/10000+1)*10000
497 endif
498
499 C Write out header info
500 write(6,*) ' Updated boundary condition data'
501 write(6,*) ' ---------------------------------'
502 write(6,*) ' Variable: ',cname
503 write(6,*) ' Description: ',cdscrip
504 write(6,*) ' Resolution: x= ',imbc,' y= ',jmbc,
505 . ' Undefined value = ',undef
506 write(6,*) ' Starting latitude = ',lat0,' longitude =',lon0
507 write(6,*) ' Data valid at these times: '
508 ndby3 = ndatebc/3
509 do n = 1, ndby3*3,3
510 write(6,1000) (nn,nymdbc(nn),nhmsbc(nn),nn=n,n+2)
511 1000 format(3(2x,i3,':',i8,2x,i8))
512 enddo
513 write(6,1000) (nn,nymdbc(nn),nhmsbc(nn),nn=ndby3*3+1,ndatebc)
514 endif
515
516 if( error ) call my_exit (101)
517
518 endif
519
520 C---------- Read SST data if necessary -------------------------------
521
522 found = .false.
523 nd = 2
524
525 c If model time is not within the times of saved sst data
526 c from previous call to getsst then read new data
527
528 timemod = float(nymdmod) + float(nhms) /1000000
529 timebc1 = float(nymdbc1) + float(nhmsbc1)/1000000
530 timebc2 = float(nymdbc2) + float(nhmsbc2)/1000000
531 if (timemod .lt. timebc1 .or. timemod .ge. timebc2) then
532
533 do while (.not.found .and. nd .le. ndatebc)
534 timebc2 = float(nymdbc(nd)) + float(nhmsbc(nd))/1000000
535 if (timebc2 .gt. timemod) then
536 nymdbc1 = nymdbc(nd-1)
537 nymdbc2 = nymdbc(nd)
538 nhmsbc1 = nhmsbc(nd-1)
539 nhmsbc2 = nhmsbc(nd)
540 call bcdata (iunit,imbc,jmbc,nd,nd+1,sstbc1,sstbc2)
541 found = .true.
542 else
543 nd = nd + 1
544 endif
545 enddo
546
547 c Otherwise the data from the last time in getsst surrounds the
548 c current model time.
549
550 else
551 found = .true.
552 endif
553
554 if (.not.found) then
555 print *, 'STOP: Could not find SST dates for model time.'
556 call my_finalize
557 call my_exit (101)
558 endif
559
560 C---------- Interpolate SST data ------------------------------------
561
562 call interp_time(nymdmod,nhms,nymdbc1,nhmsbc1,nymdbc2,nhmsbc2,
563 . fac1,fac2)
564
565 do j = jm1,jm2
566 do i = im1,im2
567 sst(i,j,bi,bj) = sstbc1(i+bislot,j+bjslot)*fac1
568 . + sstbc2(i+bislot,j+bjslot)*fac2
569 enddo
570 enddo
571
572 return
573 end
574
575 subroutine bcdata (iunit,im,jm,nrec1,nrec2,field1,field2)
576 C************************************************************************
577 C
578 C!ROUTINE: BCDATA
579 C!DESCRIPTION: BCDATA reads the data from the file assigned to the
580 C! passed unit number and returns data from the two times
581 C! surrounding the current model time. The two record
582 C! postitions are not assumed to be next to each other.
583 C
584 C!INPUT PARAMETERS:
585 C! im number of x points
586 C! im number of x points
587 C! nrec1 record number of the time before the model time
588 C! nrec2 record number of the time after the model time
589 C
590 C!OUTPUT PARAMETERS:
591 C! field1(im,jm) data field before the model time
592 C! field2(im,jm) data field after the model time
593 C
594 C--------------------------------------------------------------------------
595 implicit none
596
597 integer iunit,im,jm,nrec1,nrec2
598
599 _RL field1(im,jm)
600 _RL field2(im,jm)
601
602 integer i,j
603 real*4 f1(im,jm), f2(im,jm)
604
605 C--------- Read file -----------------------------------------------
606 read(iunit,rec=nrec1) f1
607 read(iunit,rec=nrec2) f2
608
609 #ifdef _BYTESWAPIO
610 call MDS_BYTESWAPR4( im*jm, f1)
611 call MDS_BYTESWAPR4( im*jm, f2)
612 #endif
613 do j=1,jm
614 do i=1,im
615 field1(i,j) = f1(i,j)
616 field2(i,j) = f2(i,j)
617 enddo
618 enddo
619
620 return
621 end
622 subroutine bcheader (iunit, ndmax, nrec,
623 . cname, cdscrip, im, jm, lat0, lon0, ndatebc,
624 . nymdbc, nhmsbc, undef, error)
625 C************************************************************************
626 C
627 C!ROUTINE: BCHEADER
628 C!DESCRIPTION: BCHEADER reads the header from a file and returns the info.
629 C
630 C!INPUT PARAMETERS:
631 C! iunit unit number assigned to the data file
632 C! ndmax maximum number of date/times of the data
633 C! nrec record number of the header info (or assume 1??)
634 C
635 C!OUTPUT PARAMETERS:
636 C! cname name of the data in the file header
637 C! cdscrip description of the data in the file header
638 C! im number of x points
639 C! jm number of y points
640 C! lat0 starting latitude for the data grid
641 C! lon0 starting longitude for the data grid
642 C! ndatebc number of date/times of the data in the file
643 C! nymdbc(ndmax) array of dates for the data including century
644 C! nhmsbc(ndmax) array of times for the data
645 C! undef value for undefined values in the data
646 C! error logical TRUE if dataset problem
647 C
648 C--------------------------------------------------------------------------
649 implicit none
650
651 integer iunit, ndmax, nrec
652
653 character*8 cname
654 character*80 cdscrip
655 character*112 dummy112
656 integer im,jm,ndatebc,nymdbc(ndmax),nhmsbc(ndmax)
657 _RL lat0,lon0,undef
658 logical error
659
660 integer i
661 integer*4 im_32,jm_32
662 integer*4 ndatebc_32,nhmsbc_32(ndmax),nymdbc_32(ndmax)
663 real*4 lat0_32,lon0_32,undef_32
664
665 C--------- Read file -----------------------------------------------
666
667 read(iunit,rec=nrec,err=500) cname, cdscrip,
668 . im_32, jm_32, lat0_32, lon0_32,
669 . ndatebc_32, undef_32
670
671 #ifdef _BYTESWAPIO
672 call MDS_BYTESWAPI4( 1, im_32)
673 call MDS_BYTESWAPI4( 1, jm_32)
674 call MDS_BYTESWAPR4( 1, lat0_32)
675 call MDS_BYTESWAPR4( 1, lon0_32)
676 call MDS_BYTESWAPI4( 1, ndatebc_32)
677 call MDS_BYTESWAPR4( 1, undef_32)
678 #endif
679
680 read(iunit,rec=nrec,err=500) dummy112,
681 . (nymdbc_32(i), nhmsbc_32(i), i=1,ndatebc_32)
682
683 im = im_32
684 jm = jm_32
685 lat0 = lat0_32
686 lon0 = lon0_32
687 undef = undef_32
688
689 ndatebc = ndatebc_32
690 do i=1,ndatebc
691 #ifdef _BYTESWAPIO
692 call MDS_BYTESWAPI4( 1, nymdbc_32(i))
693 call MDS_BYTESWAPI4( 1, nhmsbc_32(i))
694 #endif
695 nymdbc(i) = nymdbc_32(i)
696 nhmsbc(i) = nhmsbc_32(i)
697 enddo
698
699 return
700 500 continue
701 print *, 'Error reading boundary condition from unit ',iunit
702 error = .true.
703 return
704 end

  ViewVC Help
Powered by ViewVC 1.1.22