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 trunk/NEMO/OPA_SRC/OBC – NEMO

source: trunk/NEMO/OPA_SRC/OBC/obcdta.F90 @ 2031

Last change on this file since 2031 was 2031, checked in by rblod, 14 years ago

Fix issues with OBC and AGRIF, see ticket #688 and #692

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