/[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.26 - (show annotations) (download)
Tue Mar 13 01:38:18 2007 UTC (17 years, 3 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint58w_post, checkpoint58x_post, checkpoint59d, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59, checkpoint58y_post
Changes since 1.25: +15 -7 lines
Finish code for option for use of observed ssts

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

  ViewVC Help
Powered by ViewVC 1.1.22