/[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.23 - (show annotations) (download)
Tue Jun 14 18:14:21 2005 UTC (18 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint58e_post, checkpoint57v_post, checkpoint58u_post, checkpoint57m_post, checkpoint57s_post, checkpoint57k_post, checkpoint58r_post, checkpoint57i_post, checkpoint57y_post, checkpoint58n_post, checkpoint58t_post, checkpoint58h_post, checkpoint57y_pre, checkpoint58q_post, checkpoint58j_post, checkpoint57r_post, checkpoint58, checkpoint58f_post, checkpoint57x_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, checkpoint58o_post, checkpoint57z_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint58p_post, checkpoint57j_post, checkpoint58b_post, checkpoint58m_post, checkpoint57l_post
Changes since 1.22: +63 -40 lines
Fix to do IO only when really needed - other little stuff too

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

  ViewVC Help
Powered by ViewVC 1.1.22