New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
obcdta.F90 in tags/nemo_v3_2/nemo_v3_2/NEMO/OPA_SRC/OBC – NEMO

source: tags/nemo_v3_2/nemo_v3_2/NEMO/OPA_SRC/OBC/obcdta.F90 @ 1878

Last change on this file since 1878 was 1878, checked in by flavoni, 14 years ago

initial test for nemogcm

File size: 57.8 KB
Line 
1MODULE obcdta
2  !!==============================================================================
3  !!                            ***  MODULE obcdta  ***
4  !! Open boundary data : read the data for the open boundaries.
5  !!==============================================================================
6#if defined key_obc
7  !!------------------------------------------------------------------------------
8  !!   'key_obc'         :                                Open Boundary Conditions
9  !!------------------------------------------------------------------------------
10  !!   obc_dta           : read u, v, t, s data along each open boundary
11  !!------------------------------------------------------------------------------
12  !! * Modules used
13  USE oce             ! ocean dynamics and tracers
14  USE dom_oce         ! ocean space and time domain
15  USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
16  USE phycst          ! physical constants
17  USE obc_oce         ! ocean open boundary conditions
18  USE in_out_manager  ! I/O logical units
19  USE lib_mpp         ! distributed memory computing
20  USE dynspg_oce
21  USE ioipsl          ! now only for  ymds2ju function
22  USE iom             !
23
24  IMPLICIT NONE
25  PRIVATE
26
27  !! * Accessibility
28  PUBLIC obc_dta      ! routines called by step.F90
29  PUBLIC obc_dta_bt   ! routines called by dynspg_ts.F90
30
31  !! * Shared module variables
32  REAL(wp),  DIMENSION(2)              ::  zjcnes_obc   !
33  REAL(wp),  DIMENSION(:), ALLOCATABLE :: ztcobc
34  REAL(wp) :: rdt_obc
35  REAL(wp) :: zjcnes
36  INTEGER :: imm0, iyy0, idd0, iyy, imm, idd
37  INTEGER :: nt_a=2, nt_b=1, itobc, ndate0_cnes, nday_year0
38  INTEGER ::  itobce, itobcw, itobcs, itobcn, itobc_b  ! number of time steps in OBC files
39
40  INTEGER ::   &
41       ntobc      , &     !:  where we are in the obc file
42       ntobc_b    , &     !:  first record used
43       ntobc_a            !:  second record used
44
45  CHARACTER (len=40) :: &    ! name of data files
46       cl_obc_eTS   , cl_obc_eU,  &
47       cl_obc_wTS   , cl_obc_wU,  &
48       cl_obc_nTS   , cl_obc_nV,  &
49       cl_obc_sTS   , cl_obc_sV
50
51  ! arrays used for interpolating time dependent data on the boundaries
52  REAL(wp), DIMENSION(jpj,jpk,jptobc) :: uedta, vedta, tedta, sedta    ! East
53  REAL(wp), DIMENSION(jpj,jpk,jptobc) :: uwdta, vwdta, twdta, swdta    ! West
54  REAL(wp), DIMENSION(jpi,jpk,jptobc) :: undta, vndta, tndta, sndta    ! North
55  REAL(wp), DIMENSION(jpi,jpk,jptobc) :: usdta, vsdta, tsdta, ssdta    ! South
56
57  LOGICAL, DIMENSION (jpj,jpk ) :: ltemsk=.TRUE., luemsk=.TRUE., lvemsk=.TRUE.  ! boolean msks
58  LOGICAL, DIMENSION (jpj,jpk ) :: ltwmsk=.TRUE., luwmsk=.TRUE., lvwmsk=.TRUE.  ! used for outliers
59  LOGICAL, DIMENSION (jpi,jpk ) :: ltnmsk=.TRUE., lunmsk=.TRUE., lvnmsk=.TRUE.  ! checks
60  LOGICAL, DIMENSION (jpi,jpk ) :: ltsmsk=.TRUE., lusmsk=.TRUE., lvsmsk=.TRUE.
61
62  !! * Substitutions
63#  include "obc_vectopt_loop_substitute.h90"
64   !!----------------------------------------------------------------------
65   !!   OPA 9.0 , LOCEAN-IPSL (2006)
66   !! $Id: obcdta.F90 1732 2009-11-16 15:15:24Z smasson $
67   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
68   !!----------------------------------------------------------------------
69
70CONTAINS
71
72  SUBROUTINE obc_dta( kt )
73    !!---------------------------------------------------------------------------
74    !!                      ***  SUBROUTINE obc_dta  ***
75    !!                   
76    !! ** Purpose :   Find the climatological  boundary arrays for the specified date,
77    !!                The boundary arrays are netcdf files. Three possible cases:
78    !!                - one time frame only in the file (time dimension = 1).
79    !!                in that case the boundary data does not change in time.
80    !!                - many time frames. In that case,  if we have 12 frames
81    !!                we assume monthly fields.
82    !!                Else, we assume that time_counter is in seconds
83    !!                since the beginning of either the current year or a reference
84    !!                year given in the namelist.
85    !!                (no check is done so far but one would have to check the "unit"
86    !!                 attribute of variable time_counter).
87    !!
88    !!
89    !! History :
90    !!        !  98-05 (J.M. Molines) Original code
91    !!   8.5  !  02-10 (C. Talandier, A-M. Treguier) Free surface, F90
92    !!
93    !!   9.0  !  04-06 (F. Durand, A-M. Treguier) Netcdf BC files on input
94    !!        !  2007-2008 (C. Langlais, P. Mathiot, J.M. Molines) high frequency boundaries data
95    !!---------------------------------------------------------------------------
96    !! * Arguments
97    INTEGER, INTENT( in ) ::   kt          ! ocean time-step index
98
99    !! * Local declarations
100    INTEGER ::   ji, jj, jk, ii, ij ,it  ! dummy loop indices
101    INTEGER ::  ikprint                                ! frequency for printouts.
102    INTEGER, SAVE :: immfile, iyyfile                     !
103    INTEGER :: nt              !  record indices (incrementation)
104    INTEGER :: istop           ! local error check
105
106    REAL(wp) ::   zxy, znum, zden ! time interpolation weight
107
108    ! variables for the julian day calculation
109    INTEGER :: iyear, imonth, iday
110    REAL(wp) :: zsec , zjulian, zjuliancnes   
111
112    ! IOM STUFF
113    INTEGER ::  idvar, id_e, id_w, id_n, id_s, id_x       ! file identifiers
114    INTEGER, DIMENSION(1)  :: itmp
115    CHARACTER(LEN=25) :: cl_vname
116
117    !!---------------------------------------------------------------------------
118
119    ! 0.  initialisation :
120    ! --------------------
121    IF ( kt == nit000  )  CALL obc_dta_ini ( kt )
122    IF ( nobc_dta == 0 )  RETURN   ! already done in obc_dta_ini
123    IF ( itobc == 1    )  RETURN   ! case of only one time frame in file done in obc_dta_ini
124
125    ! in the following code, we assume that obc data are read from files, with more than 1 time frame in it
126
127    iyyfile=iyy ; immfile = 00  ! set component of the current file name
128    IF ( cffile /= 'annual') immfile = imm   !
129    IF ( ln_obc_clim       ) iyyfile = 0000  ! assume that climatological files are labeled y0000
130
131    ! 1. Synchronize time of run with time of data files
132    !---------------------------------------------------
133    ! nday_year is the day number in the current year ( 1 for 01/01 )
134    zsec=MOD( (kt-nit000)*rdt - (nday_year - nday_year0 )*rday, rday ) ! number of seconds in the current day
135    IF (ln_obc_clim)  THEN
136      zjcnes = nday_year - 1  + zsec/rday
137    ELSE
138      zjcnes = zjcnes + rdt/rday
139    ENDIF
140
141    ! look for 'before' record number in the current file
142    ntobc = nrecbef ()  ! this function return the record number for 'before', relative to zjcnes
143
144    IF (MOD(kt-1,10)==0) THEN
145       IF (lwp) WRITE(numout,*) 'kt= ',kt,' zjcnes =', zjcnes,' ndastp =',ndastp, 'mm =',imm 
146    END IF
147
148    ! 2. read a new data if necessary
149    !--------------------------------
150    IF ( ntobc /= ntobc_b ) THEN
151    ! we need to read the 'after' record
152    ! swap working index:
153    nt=nt_b ; nt_b=nt_a ; nt_a=nt
154    ntobc_b = ntobc
155
156    ! new record number :
157    ntobc_a = ntobc_a + 1 
158
159    ! all tricky things related to record number, changing files etc... are managed by obc_read
160     
161    CALL obc_read (kt, nt_a, ntobc_a, iyyfile, immfile )
162
163    ! update zjcnes_obc
164    zjcnes_obc(nt_b)= ztcobc(ntobc_b)
165    zjcnes_obc(nt_a)= ztcobc(ntobc_a)
166    ENDIF
167
168    ! 3.   interpolation at each time step
169    ! ------------------------------------
170    IF ( ln_obc_clim) THEN
171      znum= MOD(zjcnes           - zjcnes_obc(nt_b), REAL(nyear_len(1),wp) ) ; IF ( znum < 0 ) znum = znum + REAL(nyear_len(1),wp)
172      zden= MOD(zjcnes_obc(nt_a) - zjcnes_obc(nt_b), REAL(nyear_len(1),wp) ) ; IF ( zden < 0 ) zden = zden + REAL(nyear_len(1),wp)
173    ELSE
174      znum= zjcnes           - zjcnes_obc(nt_b)
175      zden= zjcnes_obc(nt_a) - zjcnes_obc(nt_b)
176    ENDIF
177    zxy = znum / zden
178
179    IF( lp_obc_east ) THEN
180       !  fills sfoe, tfoe, ufoe ,vfoe
181       sfoe(:,:) = zxy * sedta (:,:,nt_a) + (1. - zxy)*sedta(:,:,nt_b)
182       tfoe(:,:) = zxy * tedta (:,:,nt_a) + (1. - zxy)*tedta(:,:,nt_b)
183       ufoe(:,:) = zxy * uedta (:,:,nt_a) + (1. - zxy)*uedta(:,:,nt_b)
184       vfoe(:,:) = zxy * vedta (:,:,nt_a) + (1. - zxy)*vedta(:,:,nt_b)
185    ENDIF
186
187    IF( lp_obc_west) THEN
188       !  fills sfow, tfow, ufow ,vfow
189       sfow(:,:) = zxy * swdta (:,:,nt_a) + (1. - zxy)*swdta(:,:,nt_b)
190       tfow(:,:) = zxy * twdta (:,:,nt_a) + (1. - zxy)*twdta(:,:,nt_b)
191       ufow(:,:) = zxy * uwdta (:,:,nt_a) + (1. - zxy)*uwdta(:,:,nt_b)
192       vfow(:,:) = zxy * vwdta (:,:,nt_a) + (1. - zxy)*vwdta(:,:,nt_b)
193    ENDIF
194
195    IF( lp_obc_north) THEN
196       !  fills sfon, tfon, ufon ,vfon
197       sfon(:,:) = zxy * sndta (:,:,nt_a) + (1. - zxy)*sndta(:,:,nt_b)
198       tfon(:,:) = zxy * tndta (:,:,nt_a) + (1. - zxy)*tndta(:,:,nt_b)
199       ufon(:,:) = zxy * undta (:,:,nt_a) + (1. - zxy)*undta(:,:,nt_b)
200       vfon(:,:) = zxy * vndta (:,:,nt_a) + (1. - zxy)*vndta(:,:,nt_b)
201    ENDIF
202
203    IF( lp_obc_south) THEN
204       !  fills sfos, tfos, ufos ,vfos
205       sfos(:,:) = zxy * ssdta (:,:,nt_a) + (1. - zxy)*ssdta(:,:,nt_b)
206       tfos(:,:) = zxy * tsdta (:,:,nt_a) + (1. - zxy)*tsdta(:,:,nt_b)
207       ufos(:,:) = zxy * usdta (:,:,nt_a) + (1. - zxy)*usdta(:,:,nt_b)
208       vfos(:,:) = zxy * vsdta (:,:,nt_a) + (1. - zxy)*vsdta(:,:,nt_b)
209    ENDIF
210  END SUBROUTINE obc_dta
211
212
213  SUBROUTINE obc_dta_ini (kt)
214    !!-----------------------------------------------------------------------------
215    !!                       ***  SUBROUTINE obc_dta_ini  ***
216    !!
217    !! ** Purpose :
218    !!      When obc_dta first call, realize some data initialization
219    !!
220    !! ** Method :
221    !!
222    !! History :
223    !!   9.0  ! 07-10 (J.M. Molines )
224    !!----------------------------------------------------------------------------
225    !! * Argument
226    INTEGER, INTENT(in)  :: kt      ! ocean time-step index
227
228    !! * Local declarations
229    INTEGER ::   ji,jj, it   ! dummy loop indices
230
231    REAL(wp) ::   zxy                                    ! time interpolation weight
232
233    INTEGER ::  ikprint                                ! frequency for printouts.
234
235    INTEGER, SAVE :: immfile, iyyfile                     !
236    INTEGER :: nt              !  record indices (incrementation)
237    INTEGER :: istop           ! local error check
238
239    ! variables for the julian day calculation
240    INTEGER :: iyear, imonth, iday
241    REAL(wp) :: zsec , zjulian, zjuliancnes   
242
243    ! IOM STUFF
244    INTEGER ::  idvar, id_e, id_w, id_n, id_s, id_x       ! file identifiers
245    INTEGER, DIMENSION(1)  :: itmp
246    CHARACTER(LEN=25) :: cl_vname
247
248    IF(lwp) WRITE(numout,*)
249    IF(lwp) WRITE(numout,*)  'obc_dta : find boundary data'
250    IF(lwp) WRITE(numout,*)  '~~~~~~~'
251    IF (lwp) THEN
252       IF ( nobc_dta == 0 ) THEN
253          WRITE(numout,*)  '          OBC data taken from initial conditions.'
254       ELSE     
255          WRITE(numout,*)  '          OBC data taken from netcdf files.'
256       ENDIF
257    ENDIF
258    nday_year0 = nday_year  ! to remember the day when kt=nit000
259
260    sedta(:,:,:) = 0.e0 ; tedta(:,:,:) = 0.e0 ; uedta(:,:,:) = 0.e0 ; vedta(:,:,:) = 0.e0 ! East
261    swdta(:,:,:) = 0.e0 ; twdta(:,:,:) = 0.e0 ; uwdta(:,:,:) = 0.e0 ; vwdta(:,:,:) = 0.e0 ! West
262    sndta(:,:,:) = 0.e0 ; tndta(:,:,:) = 0.e0 ; undta(:,:,:) = 0.e0 ; vndta(:,:,:) = 0.e0 ! North
263    ssdta(:,:,:) = 0.e0 ; tsdta(:,:,:) = 0.e0 ; usdta(:,:,:) = 0.e0 ; vsdta(:,:,:) = 0.e0 ! South
264
265    sfoe(:,:) = 0.e0  ; tfoe(:,:) = 0.e0 ; ufoe(:,:) = 0.e0 ; vfoe(:,:) = 0.e0   ! East
266    sfow(:,:) = 0.e0  ; tfow(:,:) = 0.e0 ; ufow(:,:) = 0.e0 ; vfow(:,:) = 0.e0   ! West
267    sfon(:,:) = 0.e0  ; tfon(:,:) = 0.e0 ; ufon(:,:) = 0.e0 ; vfon(:,:) = 0.e0   ! North
268    sfos(:,:) = 0.e0  ; tfos(:,:) = 0.e0 ; ufos(:,:) = 0.e0 ; vfos(:,:) = 0.e0   ! South
269
270    IF (nobc_dta == 0 ) THEN   ! boundary data are the initial data of this run (set only at nit000)
271       IF (lp_obc_east) THEN  ! East
272          DO ji = nie0 , nie1   
273             sfoe(nje0p1:nje1m1,:) = temsk(nje0p1:nje1m1,:) * sn (ji+1 , nje0p1:nje1m1 , :)
274             tfoe(nje0p1:nje1m1,:) = temsk(nje0p1:nje1m1,:) * tn (ji+1 , nje0p1:nje1m1 , :)
275             ufoe(nje0p1:nje1m1,:) = uemsk(nje0p1:nje1m1,:) * un (ji   , nje0p1:nje1m1 , :)
276             vfoe(nje0p1:nje1m1,:) = vemsk(nje0p1:nje1m1,:) * vn (ji+1 , nje0p1:nje1m1 , :)
277          END DO
278       ENDIF
279
280       IF (lp_obc_west) THEN  ! West
281          DO ji = niw0 , niw1   
282             sfow(njw0p1:njw1m1,:) = twmsk(njw0p1:njw1m1,:) * sn (ji , njw0p1:njw1m1 , :)
283             tfow(njw0p1:njw1m1,:) = twmsk(njw0p1:njw1m1,:) * tn (ji , njw0p1:njw1m1 , :)
284             ufow(njw0p1:njw1m1,:) = uwmsk(njw0p1:njw1m1,:) * un (ji , njw0p1:njw1m1 , :)
285             vfow(njw0p1:njw1m1,:) = vwmsk(njw0p1:njw1m1,:) * vn (ji , njw0p1:njw1m1 , :)
286          END DO
287       ENDIF
288
289       IF (lp_obc_north) THEN ! North
290          DO jj = njn0 , njn1
291             sfon(nin0p1:nin1m1,:) = tnmsk(nin0p1:nin1m1,:) * sn (nin0p1:nin1m1 , jj+1 , :)
292             tfon(nin0p1:nin1m1,:) = tnmsk(nin0p1:nin1m1,:) * tn (nin0p1:nin1m1 , jj+1 , :)
293             ufon(nin0p1:nin1m1,:) = unmsk(nin0p1:nin1m1,:) * un (nin0p1:nin1m1 , jj+1 , :)
294             vfon(nin0p1:nin1m1,:) = vnmsk(nin0p1:nin1m1,:) * vn (nin0p1:nin1m1 , jj   , :)
295          END DO
296       ENDIF
297
298       IF (lp_obc_south) THEN ! South
299          DO jj = njs0 , njs1
300             sfos(nis0p1:nis1m1,:) = tsmsk(nis0p1:nis1m1,:) * sn (nis0p1:nis1m1 , jj , :)
301             tfos(nis0p1:nis1m1,:) = tsmsk(nis0p1:nis1m1,:) * tn (nis0p1:nis1m1 , jj , :)
302             ufos(nis0p1:nis1m1,:) = usmsk(nis0p1:nis1m1,:) * un (nis0p1:nis1m1 , jj , :)
303             vfos(nis0p1:nis1m1,:) = vsmsk(nis0p1:nis1m1,:) * vn (nis0p1:nis1m1 , jj , :)
304          END DO
305       ENDIF
306       RETURN  ! exit the routine all is done
307    ENDIF  ! nobc_dta = 0
308
309!!!! In the following OBC data are read from files.
310    ! all logical-mask are initialzed to true when declared
311    WHERE ( temsk == 0 ) ltemsk=.FALSE. 
312    WHERE ( uemsk == 0 ) luemsk=.FALSE. 
313    WHERE ( vemsk == 0 ) lvemsk=.FALSE. 
314
315    WHERE ( twmsk == 0 ) ltwmsk=.FALSE. 
316    WHERE ( uwmsk == 0 ) luwmsk=.FALSE. 
317    WHERE ( vwmsk == 0 ) lvwmsk=.FALSE. 
318
319    WHERE ( tnmsk == 0 ) ltnmsk=.FALSE. 
320    WHERE ( unmsk == 0 ) lunmsk=.FALSE. 
321    WHERE ( vnmsk == 0 ) lvnmsk=.FALSE. 
322
323    WHERE ( tsmsk == 0 ) ltsmsk=.FALSE. 
324    WHERE ( usmsk == 0 ) lusmsk=.FALSE. 
325    WHERE ( vsmsk == 0 ) lvsmsk=.FALSE. 
326
327    iyear=1950;  imonth=01; iday=01;  zsec=0. 
328    ! zjuliancnes : julian day corresonding  to  01/01/1950
329    CALL ymds2ju(iyear, imonth, iday,zsec , zjuliancnes)
330
331    !current year and curent month
332    iyy=INT(ndastp/10000) ; imm=INT((ndastp -iyy*10000)/100) ; idd=(ndastp-iyy*10000-imm*100)
333    IF (iyy <  1900)  iyy = iyy+1900  ! always assume that years are on 4 digits.
334    CALL ymds2ju(iyy, imm, idd ,zsec , zjulian)
335    ndate0_cnes = zjulian - zjuliancnes   ! jcnes day when call to obc_dta_ini
336
337    iyyfile=iyy ; immfile=0  ! set component of the current file name
338    IF ( cffile /= 'annual') immfile=imm
339    IF ( ln_obc_clim) iyyfile = 0  ! assume that climatological files are labeled y0000
340
341    CALL obc_dta_chktime ( iyyfile, immfile )
342
343    IF ( itobc == 1 ) THEN 
344       ! in this case we will provide boundary data only once.
345       nt_a=1 ; ntobc_a=1
346       CALL obc_read (nit000, nt_a, ntobc_a, iyyfile, immfile) 
347       IF( lp_obc_east ) THEN
348          !  fills sfoe, tfoe, ufoe ,vfoe
349          sfoe(:,:) =  sedta (:,:,1) ; tfoe(:,:) =  tedta (:,:,1)
350          ufoe(:,:) =  uedta (:,:,1) ; vfoe(:,:) =  vedta (:,:,1)
351       ENDIF
352
353       IF( lp_obc_west) THEN
354          !  fills sfow, tfow, ufow ,vfow
355          sfow(:,:) =  swdta (:,:,1) ; tfow(:,:) =  twdta (:,:,1)
356          ufow(:,:) =  uwdta (:,:,1) ; vfow(:,:) =  vwdta (:,:,1)
357       ENDIF
358
359       IF( lp_obc_north) THEN
360          !  fills sfon, tfon, ufon ,vfon
361          sfon(:,:) =  sndta (:,:,1) ; tfon(:,:) =  tndta (:,:,1)
362          ufon(:,:) =  undta (:,:,1) ; vfon(:,:) =  vndta (:,:,1)
363       ENDIF
364
365       IF( lp_obc_south) THEN
366          !  fills sfos, tfos, ufos ,vfos
367          sfos(:,:) =  ssdta (:,:,1) ; tfos(:,:) =  tsdta (:,:,1)
368          ufos(:,:) =  usdta (:,:,1) ; vfos(:,:) =  vsdta (:,:,1)
369       ENDIF
370       RETURN  ! we go out of obc_dta_ini -------------------------------------->>>>>
371    ENDIF
372
373    ! nday_year is the day number in the current year ( 1 for 01/01 )
374    ! we suppose that we always start from the begining of a day
375    !    zsec=MOD( (kt-nit000)*rdt - (nday_year - nday_year0 )*rday, rday ) ! number of seconds in the current day
376    zsec=0.e0  ! here, kt=nit000, nday_year = ndat_year0
377
378    IF (ln_obc_clim)  THEN
379      zjcnes = nday_year - 1  + zsec/rday  ! for clim file time is in days in a year
380    ELSE
381      zjcnes = ndate0_cnes + (nday_year - nday_year0 ) + zsec/rday
382    ENDIF
383
384    ! look for 'before' record number in the current file
385    ntobc = nrecbef ()
386
387    IF (lwp) WRITE(numout,*) 'obc files frequency :',cffile
388    IF (lwp) WRITE(numout,*) ' zjcnes0 =',zjcnes,' ndastp0 =',ndastp
389    IF (lwp) WRITE(numout,*) ' annee0 ',iyy,' month0 ', imm,' day0 ', idd
390    IF (lwp) WRITE(numout,*) 'first file open :',cl_obc_nTS
391
392    ! record initialisation
393    !--------------------
394    nt_b = 1 ; nt_a = 2
395
396    ntobc_a = ntobc + 1
397    ntobc_b = ntobc
398
399    CALL obc_read (kt, nt_b, ntobc_b, iyyfile, immfile)  ! read 'before' fields
400    CALL obc_read (kt, nt_a, ntobc_a, iyyfile, immfile)  ! read 'after' fields
401
402    zjcnes_obc(nt_b)= ztcobc(ntobc_b)
403    zjcnes_obc(nt_a)= ztcobc(ntobc_a)
404    !
405  END SUBROUTINE obc_dta_ini
406
407
408  SUBROUTINE obc_dta_chktime (kyyfile, kmmfile)
409   !
410   ! check the number of time steps in the files and read ztcobc
411   !
412   ! * Arguments
413   INTEGER, INTENT(in) :: kyyfile, kmmfile
414   ! * local variables
415   INTEGER :: istop       ! error control
416   INTEGER :: ji          ! dummy loop index
417
418    INTEGER ::  idvar, id_e, id_w, id_n, id_s, id_x       ! file identifiers
419    INTEGER, DIMENSION(1)  :: itmp
420    CHARACTER(LEN=25) :: cl_vname
421
422    ntobc_a = 0; itobce =0 ; itobcw = 0; itobcn = 0; itobcs = 0
423    ! build file name
424    WRITE(cl_obc_eTS ,'("obc_east_TS_y",i4.4,"m",i2.2,".nc")'  ) kyyfile,kmmfile
425    WRITE(cl_obc_wTS ,'("obc_west_TS_y",i4.4,"m",i2.2,".nc")'  ) kyyfile,kmmfile
426    WRITE(cl_obc_nTS ,'("obc_north_TS_y",i4.4,"m",i2.2,".nc")' ) kyyfile,kmmfile
427    WRITE(cl_obc_sTS ,'("obc_south_TS_y",i4.4,"m",i2.2,".nc")' ) kyyfile,kmmfile
428
429    cl_vname = 'time_counter'
430    IF ( lp_obc_east ) THEN
431       CALL iom_open ( cl_obc_eTS , id_e )
432       idvar = iom_varid( id_e, cl_vname, kdimsz = itmp ); itobce=itmp(1)
433    ENDIF
434    IF ( lp_obc_west ) THEN
435       CALL iom_open ( cl_obc_wTS , id_w )
436       idvar = iom_varid( id_w, cl_vname, kdimsz = itmp ) ; itobcw=itmp(1)
437    ENDIF
438    IF ( lp_obc_north ) THEN
439       CALL iom_open ( cl_obc_nTS , id_n )
440       idvar = iom_varid( id_n, cl_vname, kdimsz = itmp ) ; itobcn=itmp(1)
441    ENDIF
442    IF ( lp_obc_south ) THEN
443       CALL iom_open ( cl_obc_sTS , id_s )
444       idvar = iom_varid( id_s, cl_vname, kdimsz = itmp ) ; itobcs=itmp(1)
445    ENDIF
446
447    itobc = MAX( itobce, itobcw, itobcn, itobcs )
448    istop = 0
449    IF ( lp_obc_east  .AND. itobce /= itobc ) istop = istop+1 
450    IF ( lp_obc_west  .AND. itobcw /= itobc ) istop = istop+1     
451    IF ( lp_obc_north .AND. itobcn /= itobc ) istop = istop+1
452    IF ( lp_obc_south .AND. itobcs /= itobc ) istop = istop+1 
453    nstop = nstop + istop
454
455    IF ( istop /=  0 )  THEN
456       WRITE(ctmp1,*) ' east, west, north, south: ', itobce, itobcw, itobcn, itobcs
457       CALL ctl_stop( 'obcdta : all files must have the same number of time steps', ctmp1 )
458    ENDIF
459
460    IF ( itobc == 1 ) THEN
461       IF (lwp) THEN
462          WRITE(numout,*) ' obcdta found one time step only in the OBC files'
463          IF (ln_obc_clim) THEN
464             ! OK no problem
465          ELSE
466             ln_obc_clim=.true.
467             WRITE(numout,*) ' we force ln_obc_clim to T'
468          ENDIF
469       ENDIF
470    ELSE
471       IF ( ALLOCATED(ztcobc) ) DEALLOCATE ( ztcobc )
472       ALLOCATE (ztcobc(itobc))
473       DO ji=1,1   ! use a dummy loop to read ztcobc only once
474          IF ( lp_obc_east ) THEN
475             CALL iom_gettime ( id_e, ztcobc, cl_vname ) ; CALL iom_close (id_e) ; EXIT
476          ENDIF
477          IF ( lp_obc_west ) THEN
478             CALL iom_gettime ( id_w, ztcobc, cl_vname ) ; CALL iom_close (id_w) ; EXIT
479          ENDIF
480          IF ( lp_obc_north ) THEN
481             CALL iom_gettime ( id_n, ztcobc, cl_vname ) ; CALL iom_close (id_n) ; EXIT
482          ENDIF
483          IF ( lp_obc_south ) THEN
484             CALL iom_gettime ( id_s, ztcobc, cl_vname ) ; CALL iom_close (id_s) ; EXIT
485          ENDIF
486       END DO
487       rdt_obc = ztcobc(2)-ztcobc(1)  !  just an information, not used for any computation
488       IF (lwp) WRITE(numout,*) ' obcdta found', itobc,' time steps in the OBC files'
489       IF (lwp) WRITE(numout,*) ' time step of obc data :', rdt_obc,' days'           
490     ENDIF
491     zjcnes = zjcnes - rdt/rday  ! trick : zcnes is always incremented by rdt/rday in obc_dta!
492  END SUBROUTINE obc_dta_chktime
493
494
495#if defined key_dynspg_ts || defined key_dynspg_exp
496  SUBROUTINE obc_dta_bt( kt, kbt )
497    !!---------------------------------------------------------------------------
498    !!                      ***  SUBROUTINE obc_dta  ***
499    !!
500    !! ** Purpose :   time interpolation of barotropic data for time-splitting scheme
501    !!                Data at the boundary must be in m2/s
502    !!
503    !! History :
504    !!   9.0  !  05-11 (V. garnier) Original code
505    !!---------------------------------------------------------------------------
506    !! * Arguments
507    INTEGER, INTENT( in ) ::   kt          ! ocean time-step index
508    INTEGER, INTENT( in ) ::   kbt         ! barotropic ocean time-step index
509
510    !! * Local declarations
511    INTEGER ::   ji, jj, jk, ii, ij   ! dummy loop indices
512    INTEGER ::   id_e, id_w, id_n, id_s, fid  ! file identifiers
513    INTEGER ::   itimo, iman, imois, i15
514    INTEGER ::   itobcm, itobcp, itimom, itimop
515    REAL(wp) ::  zxy
516    INTEGER ::   isrel, ikt           ! number of seconds since 1/1/1992
517    INTEGER ::   iprint              ! frequency for printouts.
518
519    !!---------------------------------------------------------------------------
520
521    ! 1.   First call: check time frames available in files.
522    ! -------------------------------------------------------
523
524    IF( kt == nit000 ) THEN
525
526       ! 1.1  Barotropic tangential velocities set to zero
527       ! -------------------------------------------------
528       IF( lp_obc_east  ) vbtfoe(:) = 0.e0
529       IF( lp_obc_west  ) vbtfow(:) = 0.e0
530       IF( lp_obc_south ) ubtfos(:) = 0.e0
531       IF( lp_obc_north ) ubtfon(:) = 0.e0
532
533       ! 1.2  Sea surface height and normal barotropic velocities set to zero
534       !                               or initial conditions if nobc_dta == 0
535       ! --------------------------------------------------------------------
536
537       IF( lp_obc_east ) THEN
538          ! initialisation to zero
539          sshedta(:,:) = 0.e0
540          ubtedta(:,:) = 0.e0
541          !                                        ! ================== !
542          IF( nobc_dta == 0 )   THEN               ! initial state used !
543             !                                     ! ================== !
544             !  Fills sedta, tedta, uedta (global arrays)
545             !  Remark: this works for njzoom = 1. Should the definition of ij include njzoom?
546             DO ji = nie0, nie1
547                DO jj = nje0p1, nje1m1
548                   ij = jj -1 + njmpp
549                   sshedta(ij,1) = sshn(ji+1,jj) * tmask(ji+1,jj,1)
550                END DO
551             END DO
552          ENDIF
553       ENDIF
554
555       IF( lp_obc_west) THEN
556          ! initialisation to zero
557          sshwdta(:,:) = 0.e0
558          ubtwdta(:,:) = 0.e0
559          !                                        ! ================== !
560          IF( nobc_dta == 0 )   THEN               ! initial state used !
561             !                                     ! ================== !
562             !  Fills swdta, twdta, uwdta (global arrays)
563             !  Remark: this works for njzoom = 1. Should the definition of ij include njzoom?
564             DO ji = niw0, niw1
565                DO jj = njw0p1, njw1m1
566                   ij = jj -1 + njmpp
567                   sshwdta(ij,1) = sshn(ji,jj) * tmask(ji,jj,1)
568                END DO
569             END DO
570          ENDIF
571       ENDIF
572
573       IF( lp_obc_north) THEN
574          ! initialisation to zero
575          sshndta(:,:) = 0.e0
576          vbtndta(:,:) = 0.e0
577          !                                        ! ================== !
578          IF( nobc_dta == 0 )   THEN               ! initial state used !
579             !                                     ! ================== !
580             !  Fills sndta, tndta, vndta (global arrays)
581             !  Remark: this works for njzoom = 1. Should the definition of ij include njzoom?
582             DO jj = njn0, njn1
583                DO ji = nin0p1, nin1m1
584                   DO jk = 1, jpkm1
585                      ii = ji -1 + nimpp
586                      vbtndta(ii,1) = vbtndta(ii,1) + vndta(ii,jk,1)*fse3v(ji,jj,jk)
587                   END DO
588                   sshndta(ii,1) = sshn(ii,jj+1) * tmask(ji,jj+1,1)
589                END DO
590             END DO
591          ENDIF
592       ENDIF
593
594       IF( lp_obc_south) THEN
595          ! initialisation to zero
596          ssdta(:,:,:) = 0.e0
597          tsdta(:,:,:) = 0.e0
598          vsdta(:,:,:) = 0.e0
599          sshsdta(:,:) = 0.e0
600          vbtsdta(:,:) = 0.e0
601          !                                        ! ================== !
602          IF( nobc_dta == 0 )   THEN               ! initial state used !
603             !                                     ! ================== !
604             !  Fills ssdta, tsdta, vsdta (global arrays)
605             !  Remark: this works for njzoom = 1. Should the definition of ij include njzoom?
606             DO jj = njs0, njs1
607                DO ji = nis0p1, nis1m1
608                   DO jk = 1, jpkm1
609                      ii = ji -1 + nimpp
610                      vbtsdta(ii,1) = vbtsdta(ii,1) + vsdta(ii,jk,1)*fse3v(ji,jj,jk)
611                   END DO
612                   sshsdta(ii,1) = sshn(ji,jj) * tmask(ii,jj,1)
613                END DO
614             END DO
615          ENDIF
616       ENDIF
617
618    ENDIF        !       END IF kt == nit000
619
620    !!------------------------------------------------------------------------------------
621    ! 2.      Initialize the time we are at. Does this every time the routine is called,
622    !         excepted when nobc_dta = 0
623    !
624    IF( nobc_dta == 0) THEN
625       itimo = 1
626       zxy   = 0
627    ELSE
628       IF(ntobc == 1) THEN
629          itimo = 1
630       ELSE IF (ntobc == 12) THEN      !   BC are monthly
631          ! we assume we have climatology in that case
632          iman  = 12
633          i15   = nday / 16
634          imois = nmonth + i15 - 1
635          IF( imois == 0 )   imois = iman
636          itimo = imois
637       ELSE
638          IF(lwp) WRITE(numout,*) 'data other than constant or monthly',kt
639          iman  = ntobc
640          itimo = FLOOR( kt*rdt / tcobc(1))
641          isrel=kt*rdt
642       ENDIF
643    ENDIF
644
645    ! 2. Read two records in the file if necessary
646    ! ---------------------------------------------
647
648    IF( nobc_dta == 1 .AND. nlecto == 1 ) THEN
649
650       IF( lp_obc_east ) THEN
651          ! ... Read datafile and set sea surface height and barotropic velocity
652          ! ... initialise the sshedta, ubtedta arrays
653          sshedta(:,0) = sshedta(:,1)
654          ubtedta(:,0) = ubtedta(:,1)
655          CALL iom_open ( 'obceast_TS.nc', id_e )
656          CALL iom_get ( id_e, jpdom_unknown, 'vossurfh', sshedta(:,1), ktime=ntobc1 )
657          CALL iom_get ( id_e, jpdom_unknown, 'vossurfh', sshedta(:,2), ktime=ntobc2 )
658          IF( lk_dynspg_ts ) THEN
659             CALL iom_get (id_e, jpdom_unknown, 'vossurfh', sshedta(:,3), ktime=ntobc2+1 )
660          ENDIF
661          CALL iom_close ( id_e )
662          !
663          CALL iom_open ( 'obceast_U.nc', id_e )
664          CALL iom_get ( id_e, jpdom_unknown, 'vozoubt', ubtedta(:,1), ktime=ntobc1 )
665          CALL iom_get ( id_e, jpdom_unknown, 'vozoubt', ubtedta(:,2), ktime=ntobc2 )
666          IF( lk_dynspg_ts ) THEN
667             CALL iom_get ( id_e, jpdom_unknown, 'vozoubt', ubtedta(:,3), ktime=ntobc2+1 )
668          ENDIF
669          CALL iom_close ( id_e )
670          ! ... Usually printout is done only once at kt = nit000, unless nprint (namelist) > 1
671          IF( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
672             WRITE(numout,*)
673             WRITE(numout,*) ' Read East OBC barotropic data records ', ntobc1, ntobc2
674             iprint = (jpjef-jpjed+1)/20 +1
675             WRITE(numout,*)
676             WRITE(numout,*) ' Sea surface height record 1'
677             CALL prihre( sshedta(:,1), jpjef-jpjed+1, 1, 1, jpjef-jpjed+1, iprint, 1, 1, -3, 1., numout )
678             WRITE(numout,*)
679             WRITE(numout,*) ' Normal transport (m2/s) record 1',iprint
680             CALL prihre( ubtedta(:,1), jpjef-jpjed+1, 1, 1, jpjef-jpjed+1, iprint, 1, 1, -3, 1., numout )
681          ENDIF
682       ENDIF
683
684       IF( lp_obc_west ) THEN
685          ! ... Read datafile and set temperature, salinity and normal velocity
686          ! ... initialise the swdta, twdta, uwdta arrays
687          sshwdta(:,0) = sshwdta(:,1)
688          ubtwdta(:,0) = ubtwdta(:,1)
689          CALL iom_open ( 'obcwest_TS.nc', id_w )
690          CALL iom_get ( id_w, jpdom_unknown, 'vossurfh', sshwdta(:,1), ktime=ntobc1 )
691          CALL iom_get ( id_w, jpdom_unknown, 'vossurfh', sshwdta(:,2), ktime=ntobc2 )
692          IF( lk_dynspg_ts ) THEN
693             CALL  ( id_w, jpdom_unknown, 'vossurfh', sshwdta(:,3), ktime=ntobc2+1 )
694          ENDIF
695          CALL iom_close ( id_w )
696          !
697          CALL iom_open ( 'obcwest_U.nc', id_w )
698          CALL iom_get ( id_w, jpdom_unknown, 'vozoubt', ubtwdta(:,1), ktime=ntobc1 )
699          CALL iom_get ( id_w, jpdom_unknown, 'vozoubt', ubtwdta(:,2), ktime=ntobc2 )
700          IF( lk_dynspg_ts ) THEN
701             CALL iom_get ( id_w, jpdom_unknown, 'vozoubt', ubtwdta(:,3), ktime=ntobc2+1 )
702          ENDIF
703          CALL iom_close ( id_w )
704          ! ... Usually printout is done only once at kt = nit000, unless nprint (namelist) > 1
705          IF( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
706             WRITE(numout,*)
707             WRITE(numout,*) ' Read West OBC barotropic data records ', ntobc1, ntobc2
708             iprint = (jpjwf-jpjwd+1)/20 +1
709             WRITE(numout,*)
710             WRITE(numout,*) ' Sea surface height record 1  - printout surface level'
711             CALL prihre( sshwdta(:,1), jpjwf-jpjwd+1, 1, 1, jpjwf-jpjwd+1, iprint, 1, 1, -3, 1., numout )
712             WRITE(numout,*)
713             WRITE(numout,*) ' Normal transport (m2/s) record 1'
714             CALL prihre( ubtwdta(:,1), jpjwf-jpjwd+1, 1, 1, jpjwf-jpjwd+1, iprint, 1, 1, -3, 1., numout )
715          ENDIF
716       ENDIF
717
718       IF( lp_obc_north) THEN
719          ! ... Read datafile and set sea surface height and barotropic velocity
720          ! ... initialise the sshndta, ubtndta arrays
721          sshndta(:,0) = sshndta(:,1)
722          vbtndta(:,0) = vbtndta(:,1)
723          CALL iom_open ( 'obcnorth_TS.nc', id_n )
724          CALL iom_get (id_n, jpdom_unknown, 'vossurfh', sshndta(:,1), ktime=ntobc1 )
725          CALL iom_get (id_n, jpdom_unknown, 'vossurfh', sshndta(:,2), ktime=ntobc2 )
726          IF( lk_dynspg_ts ) THEN
727             CALL iom_get (id_n, jpdom_unknown, 'vossurfh', sshndta(:,3), ktime=ntobc2+1 )
728          ENDIF
729          CALL iom_close ( id_n )
730
731          CALL iom_open ( 'obcnorth_V.nc', id_n )
732          CALL iom_get (id_n, jpdom_unknown, 'vomevbt', vbtndta(:,1), ktime=ntobc1 )
733          CALL iom_get (id_n, jpdom_unknown, 'vomevbt', vbtndta(:,2), ktime=ntobc2 )
734          IF( lk_dynspg_ts ) THEN
735             CALL iom_get (id_n, jpdom_unknown, 'vomevbt', vbtndta(:,3), ktime=ntobc2+1 )
736          ENDIF
737          CALL iom_close ( id_n )
738
739          ! ... Usually printout is done only once at kt = nit000, unless nprint (namelist) > 1
740          IF( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
741             WRITE(numout,*)
742             WRITE(numout,*) ' Read North OBC barotropic data records ', ntobc1, ntobc2
743             iprint = (jpinf-jpind+1)/20 +1
744             WRITE(numout,*)
745             WRITE(numout,*) ' Sea surface height record 1  - printout surface level'
746             CALL prihre( sshndta(:,1), jpinf-jpind+1, 1, 1, jpinf-jpind+1, iprint, 1, 1, -3, 1., numout )
747             WRITE(numout,*)
748             WRITE(numout,*) ' Normal transport (m2/s) record 1'
749             CALL prihre( vbtndta(:,1), jpinf-jpind+1, 1, 1, jpinf-jpind+1, iprint, 1, 1, -3, 1., numout )
750          ENDIF
751       ENDIF
752
753       IF( lp_obc_south) THEN
754          ! ... Read datafile and set sea surface height and barotropic velocity
755          ! ... initialise the sshsdta, ubtsdta arrays
756          sshsdta(:,0) = sshsdta(:,1)
757          vbtsdta(:,0) = vbtsdta(:,1)
758          CALL iom_open ( 'obcsouth_TS.nc', id_s )
759          CALL iom_get ( id_s, jpdom_unknown, 'vossurfh', sshsdta(:,1), ktime=ntobc1 )
760          CALL iom_get ( id_s, jpdom_unknown, 'vossurfh', sshsdta(:,2), ktime=ntobc2 )
761          IF( lk_dynspg_ts ) THEN
762             CALL iom_get ( id_s, jpdom_unknown, 'vossurfh', sshsdta(:,3), ktime=ntobc2+1 )
763          ENDIF
764          CALL iom_close ( id_s )
765
766          CALL iom_open ( 'obcsouth_V.nc', id_s )
767          CALL iom_get ( id_s, jpdom_unknown, 'vomevbt', vbtsdta(:,1), ktime=ntobc1 )
768          CALL iom_get ( id_s, jpdom_unknown, 'vomevbt', vbtsdta(:,2), ktime=ntobc2 )
769          IF( lk_dynspg_ts ) THEN
770             CALL iom_get ( id_s, jpdom_unknown, 'vomevbt', vbtsdta(:,3), ktime=ntobc2+1 )
771          ENDIF
772          CALL iom_close ( id_s )
773
774          ! ... Usually printout is done only once at kt = nit000, unless nprint (namelist) > 1
775          IF( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
776             WRITE(numout,*)
777             WRITE(numout,*) ' Read South OBC barotropic data records ', ntobc1, ntobc2
778             iprint = (jpisf-jpisd+1)/20 +1
779             WRITE(numout,*)
780             WRITE(numout,*) ' Sea surface height record 1  - printout surface level'
781             CALL prihre( sshsdta(:,1), jpisf-jpisd+1, 1, 1, jpisf-jpisd+1, iprint, 1, 1, -3, 1., numout )
782             WRITE(numout,*)
783             WRITE(numout,*) ' Normal transport (m2/s) record 1'
784             CALL prihre( vbtsdta(:,1), jpisf-jpisd+1, 1, 1, jpisf-jpisd+1, iprint, 1, 1, -3, 1., numout )
785          ENDIF
786       ENDIF
787
788    ENDIF        !      end of the test on the condition to read or not the files
789
790    ! 3.  Call at every time step : Linear interpolation of BCs to current time step
791    ! ----------------------------------------------------------------------
792
793    IF( lk_dynspg_ts ) THEN
794       isrel = (kt-1)*rdt + kbt*(rdt/REAL(nn_baro,wp))
795
796       IF( nobc_dta == 1 ) THEN
797          isrel = (kt-1)*rdt + kbt*(rdt/REAL(nn_baro,wp))
798          itimo  = FLOOR(  kt*rdt    / (tcobc(2)-tcobc(1)) )
799          itimom = FLOOR( (kt-1)*rdt / (tcobc(2)-tcobc(1)) )
800          itimop = FLOOR( (kt+1)*rdt / (tcobc(2)-tcobc(1)) )
801          IF( itimom == itimo .AND. itimop == itimo ) THEN
802             itobcm = ntobc1
803             itobcp = ntobc2
804
805          ELSEIF ( itimom <= itimo .AND. itimop == itimo ) THEN
806             IF(  FLOOR( isrel / (tcobc(2)-tcobc(1)) ) < itimo ) THEN
807                itobcm = ntobc1-1
808                itobcp = ntobc2-1
809             ELSE
810                itobcm = ntobc1
811                itobcp = ntobc2
812             ENDIF
813
814          ELSEIF ( itimom == itimo .AND. itimop >= itimo ) THEN
815             IF(  FLOOR( isrel / (tcobc(2)-tcobc(1)) ) < itimop ) THEN
816                itobcm = ntobc1
817                itobcp = ntobc2
818             ELSE
819                itobcm = ntobc1+1
820                itobcp = ntobc2+1
821             ENDIF
822
823          ELSEIF ( itimom == itimo-1 .AND. itimop == itimo+1 ) THEN
824             IF(  FLOOR( isrel / (tcobc(2)-tcobc(1)) ) < itimo ) THEN
825                itobcm = ntobc1-1
826                itobcp = ntobc2-1
827             ELSEIF (  FLOOR( isrel / (tcobc(2)-tcobc(1)) ) < itimop ) THEN
828                itobcm = ntobc1
829                itobcp = ntobc2
830             ELSEIF (  FLOOR( isrel / (tcobc(2)-tcobc(1)) ) == itimop ) THEN
831                itobcm = ntobc1+1
832                itobcp = ntobc2+2
833             ELSE
834                IF(lwp) WRITE(numout, *) 'obc_dta_bt: You should not have seen this print! error 1?'
835             ENDIF
836          ELSE
837             IF(lwp) WRITE(numout, *) 'obc_dta_bt: You should not have seen this print! error 2?'
838          ENDIF
839
840       ENDIF
841
842    ELSE IF( lk_dynspg_exp ) THEN
843       isrel=kt*rdt
844       itobcm = ntobc1
845       itobcp = ntobc2
846    ENDIF
847
848    IF( ntobc == 1 .OR. nobc_dta == 0 ) THEN
849       zxy = 0.e0
850    ELSE IF( ntobc == 12 ) THEN
851       zxy = FLOAT( nday + 15 - 30 * i15 ) / 30.
852    ELSE
853       zxy = (tcobc(itobcm)-FLOAT(isrel)) / (tcobc(itobcm)-tcobc(itobcp))
854    ENDIF
855
856    IF( lp_obc_east ) THEN           !  fills sshfoe, ubtfoe (local to each processor)
857       DO jj = nje0p1, nje1m1
858          ij = jj -1 + njmpp
859          sshfoe(jj) = ( zxy * sshedta(ij,2) + (1.-zxy) * sshedta(ij,1) ) * temsk(jj,1)
860          ubtfoe(jj) = ( zxy * ubtedta(ij,2) + (1.-zxy) * ubtedta(ij,1) ) * uemsk(jj,1)
861       END DO
862    ENDIF
863
864    IF( lp_obc_west) THEN            !  fills sshfow, ubtfow (local to each processor)
865       DO jj = njw0p1, njw1m1
866          ij = jj -1 + njmpp
867          sshfow(jj) = ( zxy * sshwdta(ij,2) + (1.-zxy) * sshwdta(ij,1) ) * twmsk(jj,1)
868          ubtfow(jj) = ( zxy * ubtwdta(ij,2) + (1.-zxy) * ubtwdta(ij,1) ) * uwmsk(jj,1)
869       END DO
870    ENDIF
871
872    IF( lp_obc_north) THEN           !  fills sshfon, vbtfon (local to each processor)
873       DO ji = nin0p1, nin1m1
874          ii = ji -1 + nimpp
875          sshfon(ji) = ( zxy * sshndta(ii,2) + (1.-zxy) * sshndta(ii,1) ) * tnmsk(ji,1)
876          vbtfon(ji) = ( zxy * vbtndta(ii,2) + (1.-zxy) * vbtndta(ii,1) ) * vnmsk(ji,1)
877       END DO
878    ENDIF
879
880    IF( lp_obc_south) THEN           !  fills sshfos, vbtfos (local to each processor)
881       DO ji = nis0p1, nis1m1
882          ii = ji -1 + nimpp
883          sshfos(ji) = ( zxy * sshsdta(ii,2) + (1.-zxy) * sshsdta(ii,1) ) * tsmsk(ji,1)
884          vbtfos(ji) = ( zxy * vbtsdta(ii,2) + (1.-zxy) * vbtsdta(ii,1) ) * vsmsk(ji,1)
885       END DO
886    ENDIF
887
888  END SUBROUTINE obc_dta_bt
889
890#else
891  !!-----------------------------------------------------------------------------
892  !!   Default option
893  !!-----------------------------------------------------------------------------
894  SUBROUTINE obc_dta_bt ( kt, kbt )       ! Empty routine
895    !! * Arguments
896    INTEGER,INTENT(in) :: kt
897    INTEGER, INTENT( in ) ::   kbt         ! barotropic ocean time-step index
898    WRITE(*,*) 'obc_dta_bt: You should not have seen this print! error?', kt
899    WRITE(*,*) 'obc_dta_bt: You should not have seen this print! error?', kbt
900  END SUBROUTINE obc_dta_bt
901#endif
902
903
904  !!==============================================================================
905  SUBROUTINE obc_read (kt, nt_x, ntobc_x, iyy, imm)
906     !!-------------------------------------------------------------------------
907     !!                      ***  ROUTINE obc_read  ***
908     !!
909     !! ** Purpose :  Read the boundary data in files identified by iyy and imm
910     !!               According to the validated open boundaries, return the
911     !!               following arrays :
912     !!                sedta, tedta : East OBC salinity and temperature
913     !!                uedta, vedta :   "   "  u and v velocity component     
914     !!
915     !!                swdta, twdta : West OBC salinity and temperature
916     !!                uwdta, vwdta :   "   "  u and v velocity component     
917     !!
918     !!                sndta, tndta : North OBC salinity and temperature
919     !!                undta, vndta :   "   "  u and v velocity component     
920     !!
921     !!                ssdta, tsdta : South OBC salinity and temperature
922     !!                usdta, vsdta :   "   "  u and v velocity component     
923     !!
924     !! ** Method  :  These fields are read in the record ntobc_x of the files.
925     !!               The number of records is already known. If  ntobc_x is greater
926     !!               than the number of record, this routine will look for next file,
927     !!               updating the indices (case of inter-annual obcs) or loop at the
928     !!               begining in case of climatological file (ln_obc_clim = true ).
929     !! -------------------------------------------------------------------------
930     !! History:     !  2005  ( P. Mathiot, C. Langlais ) Original code
931     !!              !  2008  ( J,M, Molines ) Use IOM and cleaning
932     !!--------------------------------------------------------------------------
933
934    ! * Arguments
935    INTEGER, INTENT( in ) :: kt, nt_x
936    INTEGER, INTENT( inout ) :: ntobc_x , iyy, imm      ! yes ! inout !
937
938    ! * Local variables
939    CHARACTER (len=40) :: &    ! file names
940         cl_obc_eTS   , cl_obc_eU,  cl_obc_eV,&
941         cl_obc_wTS   , cl_obc_wU,  cl_obc_wV,&
942         cl_obc_nTS   , cl_obc_nU,  cl_obc_nV,&
943         cl_obc_sTS   , cl_obc_sU,  cl_obc_sV
944
945    INTEGER :: ikprint
946    REAL(wp) :: zmin, zmax   ! control of boundary values
947
948    !IOM stuff
949    INTEGER :: id_e, id_w, id_n, id_s, ji, jj
950    INTEGER, DIMENSION(2) :: istart, icount
951   
952    !--------------------------------------------------------------------------
953    IF ( ntobc_x > itobc ) THEN
954      IF (ln_obc_clim) THEN  ! just loop on the same file
955        ntobc_x = 1 
956      ELSE
957        ! need to change file : it is always for an 'after' data
958        IF ( cffile == 'annual' ) THEN ! go to next year file
959          iyy = iyy + 1
960        ELSE IF ( cffile =='monthly' ) THEN  ! go to next month file
961          imm = imm + 1 
962          IF ( imm == 13 ) THEN
963            imm = 1 ; iyy = iyy + 1
964          ENDIF
965        ELSE
966         ctmp1='obcread : this type of obc file is not supported :( '
967         ctmp2=TRIM(cffile)
968         CALL ctl_stop (ctmp1, ctmp2)
969         ! cffile should be either annual or monthly ...
970        ENDIF
971       ! as the file is changed, need to update itobc etc ...
972        CALL obc_dta_chktime (iyy,imm)
973        ntobc_x = nrecbef() + 1 ! remember : this case occur for an after data
974      ENDIF
975    ENDIF
976
977    IF ( lp_obc_east ) THEN 
978       ! ... Read datafile and set temperature, salinity and normal velocity
979       ! ... initialise the sedta, tedta, uedta arrays
980       WRITE(cl_obc_eTS ,'("obc_east_TS_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
981       WRITE(cl_obc_eU  ,'("obc_east_U_y"   ,i4.4,"m",i2.2,".nc")' ) iyy,imm
982       WRITE(cl_obc_eV  ,'("obc_east_V_y"   ,i4.4,"m",i2.2,".nc")' ) iyy,imm
983       ! JMM this may change depending on the obc data format ...
984       istart(:)=(/nje0+njmpp-1,1/) ; icount(:)=(/nje1-nje0 +1,jpk/)
985       IF (lwp) WRITE(numout,*) 'read data in :', TRIM(cl_obc_eTS)
986       IF (nje1 >= nje0 ) THEN
987          CALL iom_open ( cl_obc_eTS , id_e )
988          CALL iom_get ( id_e, jpdom_unknown, 'votemper', tedta(nje0:nje1,:,nt_x), &
989               &               ktime=ntobc_x , kstart=istart, kcount= icount )
990          CALL iom_get ( id_e, jpdom_unknown, 'vosaline', sedta(nje0:nje1,:,nt_x), &
991               &               ktime=ntobc_x , kstart=istart, kcount= icount )
992          CALL iom_close (id_e)
993          !
994          CALL iom_open ( cl_obc_eU , id_e )
995          CALL iom_get  ( id_e, jpdom_unknown, 'vozocrtx', uedta(nje0:nje1,:,nt_x), &
996               &               ktime=ntobc_x , kstart=istart, kcount= icount )
997          CALL iom_close ( id_e )
998          !
999          CALL iom_open ( cl_obc_eV , id_e )
1000          CALL iom_get ( id_e, jpdom_unknown, 'vomecrty', vedta(nje0:nje1,:,nt_x), &
1001               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1002          CALL iom_close ( id_e )
1003
1004          ! mask the boundary values
1005          tedta(:,:,nt_x) = tedta(:,:,nt_x)*temsk(:,:) ;  sedta(:,:,nt_x) = sedta(:,:,nt_x)*temsk(:,:)
1006          uedta(:,:,nt_x) = uedta(:,:,nt_x)*uemsk(:,:) ;  vedta(:,:,nt_x) = vedta(:,:,nt_x)*vemsk(:,:)
1007
1008          ! check any outliers
1009          zmin=MINVAL( sedta(:,:,nt_x), mask=ltemsk ) ; zmax=MAXVAL(sedta(:,:,nt_x), mask=ltemsk)
1010          IF (  zmin < 5 .OR. zmax > 50)   THEN
1011             CALL ctl_stop('Error in sedta',' routine obcdta')
1012          ENDIF
1013          zmin=MINVAL( tedta(:,:,nt_x), mask=ltemsk ) ; zmax=MAXVAL(tedta(:,:,nt_x), mask=ltemsk)
1014          IF (  zmin < -10. .OR. zmax > 40)   THEN
1015             CALL ctl_stop('Error in tedta',' routine obcdta')
1016          ENDIF
1017          zmin=MINVAL( uedta(:,:,nt_x), mask=luemsk ) ; zmax=MAXVAL(uedta(:,:,nt_x), mask=luemsk)
1018          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1019             CALL ctl_stop('Error in uedta',' routine obcdta')
1020          ENDIF
1021          zmin=MINVAL( vedta(:,:,nt_x), mask=lvemsk ) ; zmax=MAXVAL(vedta(:,:,nt_x), mask=lvemsk)
1022          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1023             CALL ctl_stop('Error in vedta',' routine obcdta')
1024          ENDIF
1025
1026          !               Usually printout is done only once at kt = nit000, unless nprint (namelist) > 1     
1027          IF ( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
1028             WRITE(numout,*)
1029             WRITE(numout,*) ' Read East OBC data records ', ntobc_x
1030             ikprint = jpj/20 +1
1031             WRITE(numout,*) ' Temperature  record 1 - printout every 3 level'
1032             CALL prihre( tedta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1033             WRITE(numout,*)
1034             WRITE(numout,*) ' Salinity  record 1 - printout every 3 level'
1035             CALL prihre( sedta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1036             WRITE(numout,*)
1037             WRITE(numout,*) ' Normal velocity U  record 1  - printout every 3 level'
1038             CALL prihre( uedta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1039             WRITE(numout,*)
1040             WRITE(numout,*) ' Tangential velocity V  record 1  - printout every 3 level'
1041             CALL prihre( vedta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1042          ENDIF
1043       ENDIF
1044    ENDIF
1045!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1046    IF ( lp_obc_west ) THEN
1047       ! ... Read datafile and set temperature, salinity and normal velocity
1048       ! ... initialise the swdta, twdta, uwdta arrays
1049       WRITE(cl_obc_wTS ,'("obc_west_TS_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1050       WRITE(cl_obc_wU  ,'("obc_west_U_y"   ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1051       WRITE(cl_obc_wV  ,'("obc_west_V_y"   ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1052       istart(:)=(/njw0+njmpp-1,1/) ; icount(:)=(/njw1-njw0 +1,jpk/)
1053       IF (lwp) WRITE(numout,*) 'read data in :', TRIM(cl_obc_wTS)
1054
1055       IF ( njw1 >= njw0 ) THEN
1056          CALL iom_open ( cl_obc_wTS , id_w )
1057          CALL iom_get ( id_w, jpdom_unknown, 'votemper', twdta(njw0:njw1,:,nt_x), & 
1058               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1059
1060          CALL iom_get ( id_w, jpdom_unknown, 'vosaline', swdta(njw0:njw1,:,nt_x), &
1061               &               ktime=ntobc_x , kstart=istart, kcount= icount)
1062          CALL iom_close (id_w)
1063          !
1064          CALL iom_open ( cl_obc_wU , id_w )
1065          CALL iom_get  ( id_w, jpdom_unknown, 'vozocrtx', uwdta(njw0:njw1,:,nt_x),&
1066               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1067          CALL iom_close ( id_w )
1068          !
1069          CALL iom_open ( cl_obc_wV , id_w )
1070          CALL iom_get ( id_w, jpdom_unknown, 'vomecrty', vwdta(njw0:njw1,:,nt_x), &
1071               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1072          CALL iom_close ( id_w )
1073
1074          ! mask the boundary values
1075          twdta(:,:,nt_x) = twdta(:,:,nt_x)*twmsk(:,:) ;  swdta(:,:,nt_x) = swdta(:,:,nt_x)*twmsk(:,:)
1076          uwdta(:,:,nt_x) = uwdta(:,:,nt_x)*uwmsk(:,:) ;  vwdta(:,:,nt_x) = vwdta(:,:,nt_x)*vwmsk(:,:)
1077
1078          ! check any outliers
1079          zmin=MINVAL( swdta(:,:,nt_x), mask=ltwmsk ) ; zmax=MAXVAL(swdta(:,:,nt_x), mask=ltwmsk)
1080          IF (  zmin < 5 .OR. zmax > 50)   THEN
1081             CALL ctl_stop('Error in swdta',' routine obcdta')
1082          ENDIF
1083          zmin=MINVAL( twdta(:,:,nt_x), mask=ltwmsk ) ; zmax=MAXVAL(twdta(:,:,nt_x), mask=ltwmsk)
1084          IF (  zmin < -10. .OR. zmax > 40)   THEN
1085             CALL ctl_stop('Error in twdta',' routine obcdta')
1086          ENDIF
1087          zmin=MINVAL( uwdta(:,:,nt_x), mask=luwmsk ) ; zmax=MAXVAL(uwdta(:,:,nt_x), mask=luwmsk)
1088          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1089             CALL ctl_stop('Error in uwdta',' routine obcdta')
1090          ENDIF
1091          zmin=MINVAL( vwdta(:,:,nt_x), mask=lvwmsk ) ; zmax=MAXVAL(vwdta(:,:,nt_x), mask=lvwmsk)
1092          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1093             CALL ctl_stop('Error in vwdta',' routine obcdta')
1094          ENDIF
1095
1096
1097          IF ( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
1098             WRITE(numout,*)
1099             WRITE(numout,*) ' Read West OBC data records ', ntobc_x
1100             ikprint = jpj/20 +1
1101             WRITE(numout,*) ' Temperature  record 1 - printout every 3 level'
1102             CALL prihre( twdta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1103             WRITE(numout,*)
1104             WRITE(numout,*) ' Salinity  record 1 - printout every 3 level'
1105             CALL prihre( swdta(:,:,nt_x),jpj,jpk, 1, jpj, ikprint,   jpk, 1, -3, 1., numout )
1106             WRITE(numout,*)
1107             WRITE(numout,*) ' Normal velocity U  record 1  - printout every 3 level'
1108             CALL prihre( uwdta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1109             WRITE(numout,*)
1110             WRITE(numout,*) ' Tangential velocity V  record 1  - printout every 3 level'
1111             CALL prihre( vwdta(:,:,nt_x), jpj, jpk, 1, jpj, ikprint, jpk, 1, -3, 1., numout )
1112          ENDIF
1113       END IF
1114    ENDIF
1115!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1116    IF( lp_obc_north) THEN
1117       WRITE(cl_obc_nTS ,'("obc_north_TS_y" ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1118       WRITE(cl_obc_nV  ,'("obc_north_V_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1119       WRITE(cl_obc_nU  ,'("obc_north_U_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1120       istart(:)=(/nin0+nimpp-1,1/) ; icount(:)=(/nin1-nin0 +1,jpk/)
1121       IF (lwp) WRITE(numout,*) 'read data in :', TRIM(cl_obc_nTS)
1122       IF ( nin1 >= nin0 ) THEN
1123          CALL iom_open ( cl_obc_nTS , id_n )
1124          CALL iom_get ( id_n, jpdom_unknown, 'votemper', tndta(nin0:nin1,:,nt_x), &
1125               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1126          CALL iom_get ( id_n, jpdom_unknown, 'vosaline', sndta(nin0:nin1,:,nt_x), &
1127               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1128          CALL iom_close (id_n)
1129          !
1130          CALL iom_open ( cl_obc_nU , id_n )
1131          CALL iom_get  ( id_n, jpdom_unknown, 'vozocrtx', undta(nin0:nin1,:,nt_x), &
1132               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1133          CALL iom_close ( id_n )
1134          !
1135          CALL iom_open ( cl_obc_nV , id_n )
1136          CALL iom_get  ( id_n, jpdom_unknown, 'vomecrty', vndta(nin0:nin1,:,nt_x), &
1137               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1138          CALL iom_close ( id_n )
1139
1140          ! mask the boundary values
1141          tndta(:,:,nt_x) = tndta(:,:,nt_x)*tnmsk(:,:) ;  sndta(:,:,nt_x) = sndta(:,:,nt_x)*tnmsk(:,:)
1142          undta(:,:,nt_x) = undta(:,:,nt_x)*unmsk(:,:) ;  vndta(:,:,nt_x) = vndta(:,:,nt_x)*vnmsk(:,:)
1143
1144          ! check any outliers
1145          zmin=MINVAL( sndta(:,:,nt_x), mask=ltnmsk ) ; zmax=MAXVAL(sndta(:,:,nt_x), mask=ltnmsk)
1146          IF (  zmin < 5 .OR. zmax > 50)   THEN
1147             CALL ctl_stop('Error in sndta',' routine obcdta')
1148          ENDIF
1149          zmin=MINVAL( tndta(:,:,nt_x), mask=ltnmsk ) ; zmax=MAXVAL(tndta(:,:,nt_x), mask=ltnmsk)
1150          IF (  zmin < -10. .OR. zmax > 40)   THEN
1151             CALL ctl_stop('Error in tndta',' routine obcdta')
1152          ENDIF
1153          zmin=MINVAL( undta(:,:,nt_x), mask=lunmsk ) ; zmax=MAXVAL(undta(:,:,nt_x), mask=lunmsk)
1154          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1155             CALL ctl_stop('Error in undta',' routine obcdta')
1156          ENDIF
1157          zmin=MINVAL( vndta(:,:,nt_x), mask=lvnmsk ) ; zmax=MAXVAL(vndta(:,:,nt_x), mask=lvnmsk)
1158          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1159             CALL ctl_stop('Error in vndta',' routine obcdta')
1160          ENDIF
1161
1162          IF ( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
1163             WRITE(numout,*)
1164             WRITE(numout,*) ' Read North OBC data records ', ntobc_x
1165             ikprint = jpi/20 +1
1166             WRITE(numout,*) ' Temperature  record 1 - printout every 3 level'
1167             CALL prihre( tndta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1168             WRITE(numout,*)
1169             WRITE(numout,*) ' Salinity  record 1 - printout every 3 level'
1170             CALL prihre( sndta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1171             WRITE(numout,*)
1172             WRITE(numout,*) ' Normal velocity V  record 1  - printout every 3 level'
1173             CALL prihre( vndta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1174             WRITE(numout,*)
1175             WRITE(numout,*) ' Tangential  velocity U  record 1  - printout every 3 level'
1176             CALL prihre( undta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1177          ENDIF
1178       ENDIF
1179    ENDIF
1180!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1181    IF( lp_obc_south) THEN
1182       WRITE(cl_obc_sTS ,'("obc_south_TS_y" ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1183       WRITE(cl_obc_sV  ,'("obc_south_V_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1184       WRITE(cl_obc_sU  ,'("obc_south_U_y"  ,i4.4,"m",i2.2,".nc")' ) iyy,imm
1185       istart(:)=(/nis0+nimpp-1,1/) ; icount(:)=(/nis1-nis0 +1,jpk/)
1186       IF (lwp) WRITE(numout,*) 'read data in :', TRIM(cl_obc_sTS)
1187       IF ( nis1 >= nis0 ) THEN
1188          CALL iom_open ( cl_obc_sTS , id_s )
1189          CALL iom_get ( id_s, jpdom_unknown, 'votemper', tsdta(nis0:nis1,:,nt_x), &
1190               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1191          CALL iom_get ( id_s, jpdom_unknown, 'vosaline', ssdta(nis0:nis1,:,nt_x), &
1192               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1193          CALL iom_close (id_s)
1194          !
1195          CALL iom_open ( cl_obc_sU , id_s )
1196          CALL iom_get  ( id_s, jpdom_unknown, 'vozocrtx', usdta(nis0:nis1,:,nt_x), &
1197               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1198          CALL iom_close ( id_s )
1199          !
1200          CALL iom_open ( cl_obc_sV , id_s )
1201          CALL iom_get  ( id_s, jpdom_unknown, 'vomecrty', vsdta(nis0:nis1,:,nt_x), &
1202               &               ktime=ntobc_x , kstart=istart, kcount= icount )
1203          CALL iom_close ( id_s )
1204
1205          ! mask the boundary values
1206          tsdta(:,:,nt_x) = tsdta(:,:,nt_x)*tsmsk(:,:) ;  ssdta(:,:,nt_x) = ssdta(:,:,nt_x)*tsmsk(:,:)
1207          usdta(:,:,nt_x) = usdta(:,:,nt_x)*usmsk(:,:) ;  vsdta(:,:,nt_x) = vsdta(:,:,nt_x)*vsmsk(:,:)
1208
1209          ! check any outliers
1210          zmin=MINVAL( ssdta(:,:,nt_x), mask=ltsmsk ) ; zmax=MAXVAL(ssdta(:,:,nt_x), mask=ltsmsk)
1211          IF (  zmin < 5 .OR. zmax > 50)   THEN
1212             CALL ctl_stop('Error in ssdta',' routine obcdta')
1213          ENDIF
1214          zmin=MINVAL( tsdta(:,:,nt_x), mask=ltsmsk ) ; zmax=MAXVAL(tsdta(:,:,nt_x), mask=ltsmsk)
1215          IF (  zmin < -10. .OR. zmax > 40)   THEN
1216             CALL ctl_stop('Error in tsdta',' routine obcdta')
1217          ENDIF
1218          zmin=MINVAL( usdta(:,:,nt_x), mask=lusmsk ) ; zmax=MAXVAL(usdta(:,:,nt_x), mask=lusmsk)
1219          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1220             CALL ctl_stop('Error in usdta',' routine obcdta')
1221          ENDIF
1222          zmin=MINVAL( vsdta(:,:,nt_x), mask=lvsmsk ) ; zmax=MAXVAL(vsdta(:,:,nt_x), mask=lvsmsk)
1223          IF (  zmin < -5. .OR. zmax > 5.)   THEN
1224             CALL ctl_stop('Error in vsdta',' routine obcdta')
1225          ENDIF
1226
1227          IF ( lwp .AND.  ( kt == nit000 .OR. nprint /= 0 )  ) THEN
1228             WRITE(numout,*)
1229             WRITE(numout,*) ' Read South OBC data records ', ntobc_x
1230             ikprint = jpi/20 +1
1231             WRITE(numout,*) ' Temperature  record 1 - printout every 3 level'
1232             CALL prihre( tsdta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1233             WRITE(numout,*)
1234             WRITE(numout,*) ' Salinity  record 1 - printout every 3 level'
1235             CALL prihre( ssdta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1236             WRITE(numout,*)
1237             WRITE(numout,*) ' Normal velocity V  record 1  - printout every 3 level'
1238             CALL prihre( vsdta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1239             WRITE(numout,*)
1240             WRITE(numout,*) ' Tangential velocity U  record 1  - printout every 3 level'
1241             CALL prihre( usdta(:,:,nt_x), jpi, jpk, 1, jpi, ikprint, jpk, 1, -3, 1., numout )
1242          ENDIF
1243       ENDIF
1244    ENDIF
1245!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1246  END SUBROUTINE obc_read 
1247
1248  INTEGER FUNCTION nrecbef()
1249      !!-----------------------------------------------------------------------
1250      !!                     ***    FUNCTION nrecbef   ***
1251      !!
1252      !!  Purpose : - provide the before record number in files, with respect to zjcnes
1253      !!
1254      !!    History : 2008-04 : ( J.M. Molines ) Original code
1255      !!-----------------------------------------------------------------------
1256
1257      INTEGER :: it , idum
1258
1259    idum = itobc
1260    DO it =1, itobc
1261       IF ( ztcobc(it) > zjcnes ) THEN ;  idum = it - 1 ; EXIT ;  ENDIF
1262    ENDDO
1263    ! idum can be 0 (climato, before first record)
1264    IF ( idum == 0 ) THEN
1265       IF ( ln_obc_clim ) THEN
1266         idum = itobc
1267       ELSE
1268         ctmp1='obc_dta: find ntobc == 0 for  non climatological file '
1269         ctmp2='consider adding a first record in your data file '
1270         CALL ctl_stop(ctmp1, ctmp2)
1271       ENDIF
1272    ENDIF
1273    ! idum can be itobc ( zjcnes > ztcobc (itobc) )
1274    !  This is not a problem ...
1275    nrecbef = idum
1276
1277  END FUNCTION nrecbef
1278
1279#else
1280  !!------------------------------------------------------------------------------
1281  !!   default option:           Dummy module          NO Open Boundary Conditions
1282  !!------------------------------------------------------------------------------
1283CONTAINS
1284  SUBROUTINE obc_dta( kt )             ! Dummy routine
1285    INTEGER, INTENT (in) :: kt
1286    WRITE(*,*) 'obc_dta: You should not have seen this print! error?', kt
1287  END SUBROUTINE obc_dta
1288#endif
1289END MODULE obcdta
Note: See TracBrowser for help on using the repository browser.