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

source: tags/nemo_v3_2_beta/NEMO/OPA_SRC/OBC/obcdta.F90 @ 4839

Last change on this file since 4839 was 1528, checked in by rblod, 15 years ago

Suppress rigid-lid option, see ticket #486

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