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 branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBC – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBC/obcdta.F90 @ 4400

Last change on this file since 4400 was 3432, checked in by trackstand2, 12 years ago

Merge branch 'ksection_partition'

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