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.
bdydta.F90 in NEMO/branches/NERC/NEMO_4.0.4_CO9_package_tides/src/OCE/BDY – NEMO

source: NEMO/branches/NERC/NEMO_4.0.4_CO9_package_tides/src/OCE/BDY/bdydta.F90 @ 14364

Last change on this file since 14364 was 14364, checked in by deazer, 3 years ago

Just a small change as nblen not declared

File size: 42.7 KB
Line 
1MODULE bdydta
2   !!======================================================================
3   !!                       ***  MODULE bdydta  ***
4   !! Open boundary data : read the data for the unstructured open boundaries.
5   !!======================================================================
6   !! History :  1.0  !  2005-01  (J. Chanut, A. Sellar)  Original code
7   !!             -   !  2007-01  (D. Storkey) Update to use IOM module
8   !!             -   !  2007-07  (D. Storkey) add bdy_dta_fla
9   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
10   !!            3.3  !  2010-09  (E.O'Dea) modifications for Shelf configurations
11   !!            3.3  !  2010-09  (D.Storkey) add ice boundary conditions
12   !!            3.4  !  2011     (D. Storkey) rewrite in preparation for OBC-BDY merge
13   !!            3.6  !  2012-01  (C. Rousset) add ice boundary conditions for sea ice
14   !!            4.0  !  2018     (C. Rousset) SI3 compatibility
15   !!----------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------
18   !!    bdy_dta      : read external data along open boundaries from file
19   !!    bdy_dta_init : initialise arrays etc for reading of external data
20   !!----------------------------------------------------------------------
21   USE oce            ! ocean dynamics and tracers
22   USE dom_oce        ! ocean space and time domain
23   USE phycst         ! physical constants
24   USE sbcapr         ! atmospheric pressure forcing
25   USE sbctide        ! Tidal forcing or not
26   USE bdy_oce        ! ocean open boundary conditions 
27   USE bdytides       ! tidal forcing at boundaries
28#if defined key_si3
29   USE ice            ! sea-ice variables
30   USE icevar         ! redistribute ice input into categories
31#endif
32   !
33   USE lib_mpp, ONLY: ctl_stop, ctl_nam
34   USE fldread        ! read input fields
35   USE iom            ! IOM library
36   USE in_out_manager ! I/O logical units
37   USE timing         ! Timing
38
39   IMPLICIT NONE
40   PRIVATE
41
42   PUBLIC   bdy_dta          ! routine called by step.F90 and dynspg_ts.F90
43   PUBLIC   bdy_dta_init     ! routine called by nemogcm.F90
44
45   INTEGER , PARAMETER ::   jpbdyfld  = 17    ! maximum number of files to read
46   INTEGER , PARAMETER ::   jp_bdyssh = 1     !
47   INTEGER , PARAMETER ::   jp_bdyu2d = 2     !
48   INTEGER , PARAMETER ::   jp_bdyv2d = 3     !
49   INTEGER , PARAMETER ::   jp_bdyu3d = 4     !
50   INTEGER , PARAMETER ::   jp_bdyv3d = 5     !
51   INTEGER , PARAMETER ::   jp_bdytem = 6     !
52   INTEGER , PARAMETER ::   jp_bdysal = 7     !
53   INTEGER , PARAMETER ::   jp_bdya_i = 8     !
54   INTEGER , PARAMETER ::   jp_bdyh_i = 9     !
55   INTEGER , PARAMETER ::   jp_bdyh_s = 10    !
56   INTEGER , PARAMETER ::   jp_bdyt_i = 11    !
57   INTEGER , PARAMETER ::   jp_bdyt_s = 12    !
58   INTEGER , PARAMETER ::   jp_bdytsu = 13    !
59   INTEGER , PARAMETER ::   jp_bdys_i = 14    !
60   INTEGER , PARAMETER ::   jp_bdyaip = 15    !
61   INTEGER , PARAMETER ::   jp_bdyhip = 16    !
62   INTEGER , PARAMETER ::   jp_bdyhil = 17    !
63#if ! defined key_si3
64   INTEGER , PARAMETER ::   jpl = 1
65#endif
66
67!$AGRIF_DO_NOT_TREAT
68   TYPE(FLD), PUBLIC, ALLOCATABLE, DIMENSION(:,:), TARGET ::   bf   ! structure of input fields (file informations, fields read)
69!$AGRIF_END_DO_NOT_TREAT
70
71   !!----------------------------------------------------------------------
72   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
73   !! $Id$
74   !! Software governed by the CeCILL license (see ./LICENSE)
75   !!----------------------------------------------------------------------
76CONTAINS
77
78   SUBROUTINE bdy_dta( kt, kit, kt_offset )
79      !!----------------------------------------------------------------------
80      !!                   ***  SUBROUTINE bdy_dta  ***
81      !!                   
82      !! ** Purpose :   Update external data for open boundary conditions
83      !!
84      !! ** Method  :   Use fldread.F90
85      !!               
86      !!----------------------------------------------------------------------
87      INTEGER, INTENT(in)           ::   kt           ! ocean time-step index
88      INTEGER, INTENT(in), OPTIONAL ::   kit          ! subcycle time-step index (for timesplitting option)
89      INTEGER, INTENT(in), OPTIONAL ::   kt_offset    ! time offset in units of timesteps. NB. if kit
90      !                                               ! is present then units = subcycle timesteps.
91      !                                               ! kt_offset = 0 => get data at "now" time level
92      !                                               ! kt_offset = -1 => get data at "before" time level
93      !                                               ! kt_offset = +1 => get data at "after" time level
94      !                                               ! etc.
95      !
96      INTEGER ::  jbdy, jfld, jstart, jend, ib, jl    ! dummy loop indices
97      INTEGER ::  ii, ij, ik, igrd, ipl               ! local integers
98      TYPE(OBC_DATA)         , POINTER ::   dta_alias        ! short cut
99      TYPE(FLD), DIMENSION(:), POINTER ::   bf_alias
100      !!---------------------------------------------------------------------------
101      !
102      IF( ln_timing )   CALL timing_start('bdy_dta')
103      !
104      ! Initialise data arrays once for all from initial conditions where required
105      !---------------------------------------------------------------------------
106      IF( kt == nit000 .AND. .NOT.PRESENT(kit) ) THEN
107
108         ! Calculate depth-mean currents
109         !-----------------------------
110
111         DO jbdy = 1, nb_bdy
112            !
113            IF( nn_dyn2d_dta(jbdy) == 0 ) THEN
114               IF( dta_bdy(jbdy)%lneed_ssh ) THEN
115                  igrd = 1
116                  DO ib = 1, idx_bdy(jbdy)%nblenrim(igrd)   ! ssh is allocated and used only on the rim
117                     ii = idx_bdy(jbdy)%nbi(ib,igrd)
118                     ij = idx_bdy(jbdy)%nbj(ib,igrd)
119                     dta_bdy(jbdy)%ssh(ib) = sshn(ii,ij) * tmask(ii,ij,1)         
120                  END DO
121               ENDIF
122               IF( ASSOCIATED(dta_bdy(jbdy)%u2d) ) THEN   ! no SIZE with a unassociated pointer. v2d and u2d can differ on subdomain
123                  igrd = 2
124                  DO ib = 1, SIZE(dta_bdy(jbdy)%u2d)      ! u2d is used either over the whole bdy or only on the rim
125                     ii = idx_bdy(jbdy)%nbi(ib,igrd)
126                     ij = idx_bdy(jbdy)%nbj(ib,igrd)
127                     dta_bdy(jbdy)%u2d(ib) = un_b(ii,ij) * umask(ii,ij,1)         
128                  END DO
129               ENDIF
130               IF( ASSOCIATED(dta_bdy(jbdy)%v2d) ) THEN   ! no SIZE with a unassociated pointer. v2d and u2d can differ on subdomain
131                  igrd = 3
132                  DO ib = 1, SIZE(dta_bdy(jbdy)%v2d)      ! v2d is used either over the whole bdy or only on the rim
133                     ii = idx_bdy(jbdy)%nbi(ib,igrd)
134                     ij = idx_bdy(jbdy)%nbj(ib,igrd)
135                     dta_bdy(jbdy)%v2d(ib) = vn_b(ii,ij) * vmask(ii,ij,1)         
136                  END DO
137               ENDIF
138            ENDIF
139            !
140            IF( nn_dyn3d_dta(jbdy) == 0 ) THEN
141               IF( dta_bdy(jbdy)%lneed_dyn3d ) THEN
142                  igrd = 2 
143                  DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
144                     DO ik = 1, jpkm1
145                        ii = idx_bdy(jbdy)%nbi(ib,igrd)
146                        ij = idx_bdy(jbdy)%nbj(ib,igrd)
147                        dta_bdy(jbdy)%u3d(ib,ik) =  ( un(ii,ij,ik) - un_b(ii,ij) ) * umask(ii,ij,ik)         
148                     END DO
149                  END DO
150                  igrd = 3 
151                  DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
152                     DO ik = 1, jpkm1
153                        ii = idx_bdy(jbdy)%nbi(ib,igrd)
154                        ij = idx_bdy(jbdy)%nbj(ib,igrd)
155                        dta_bdy(jbdy)%v3d(ib,ik) =  ( vn(ii,ij,ik) - vn_b(ii,ij) ) * vmask(ii,ij,ik)         
156                     END DO
157                  END DO
158               ENDIF
159            ENDIF
160
161            IF( nn_tra_dta(jbdy) == 0 ) THEN
162               IF( dta_bdy(jbdy)%lneed_tra ) THEN
163                  igrd = 1 
164                  DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
165                     DO ik = 1, jpkm1
166                        ii = idx_bdy(jbdy)%nbi(ib,igrd)
167                        ij = idx_bdy(jbdy)%nbj(ib,igrd)
168                        dta_bdy(jbdy)%tem(ib,ik) = tsn(ii,ij,ik,jp_tem) * tmask(ii,ij,ik)         
169                        dta_bdy(jbdy)%sal(ib,ik) = tsn(ii,ij,ik,jp_sal) * tmask(ii,ij,ik)         
170                     END DO
171                  END DO
172               ENDIF
173            ENDIF
174
175#if defined key_si3
176            IF( nn_ice_dta(jbdy) == 0 ) THEN    ! set ice to initial values
177               IF( dta_bdy(jbdy)%lneed_ice ) THEN
178                  igrd = 1   
179                  DO jl = 1, jpl
180                     DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
181                        ii = idx_bdy(jbdy)%nbi(ib,igrd)
182                        ij = idx_bdy(jbdy)%nbj(ib,igrd)
183                        dta_bdy(jbdy)%a_i(ib,jl) =  a_i (ii,ij,jl) * tmask(ii,ij,1) 
184                        dta_bdy(jbdy)%h_i(ib,jl) =  h_i (ii,ij,jl) * tmask(ii,ij,1) 
185                        dta_bdy(jbdy)%h_s(ib,jl) =  h_s (ii,ij,jl) * tmask(ii,ij,1) 
186                        dta_bdy(jbdy)%t_i(ib,jl) =  SUM(t_i (ii,ij,:,jl)) * r1_nlay_i * tmask(ii,ij,1) 
187                        dta_bdy(jbdy)%t_s(ib,jl) =  SUM(t_s (ii,ij,:,jl)) * r1_nlay_s * tmask(ii,ij,1)
188                        dta_bdy(jbdy)%tsu(ib,jl) =  t_su(ii,ij,jl) * tmask(ii,ij,1) 
189                        dta_bdy(jbdy)%s_i(ib,jl) =  s_i (ii,ij,jl) * tmask(ii,ij,1)
190                        ! melt ponds
191                        dta_bdy(jbdy)%aip(ib,jl) =  a_ip(ii,ij,jl) * tmask(ii,ij,1) 
192                        dta_bdy(jbdy)%hip(ib,jl) =  h_ip(ii,ij,jl) * tmask(ii,ij,1) 
193                        dta_bdy(jbdy)%hil(ib,jl) =  h_il(ii,ij,jl) * tmask(ii,ij,1) 
194                     END DO
195                  END DO
196               ENDIF
197            ENDIF
198#endif
199         END DO ! jbdy
200         !
201      ENDIF ! kt == nit000
202
203      ! update external data from files
204      !--------------------------------
205
206      DO jbdy = 1, nb_bdy
207
208         dta_alias => dta_bdy(jbdy)
209         bf_alias  => bf(:,jbdy)
210
211         ! read/update all bdy data
212         ! ------------------------
213         CALL fld_read( kt, 1, bf_alias, kit = kit, kt_offset = kt_offset )
214
215         ! apply some corrections in some specific cases...
216         ! --------------------------------------------------
217         !
218         ! if runoff condition: change river flow we read (in m3/s) into barotropic velocity (m/s)
219         IF( cn_tra(jbdy) == 'runoff' ) THEN   ! runoff
220            !
221            IF( ASSOCIATED(dta_bdy(jbdy)%u2d) ) THEN   ! no SIZE with a unassociated pointer. v2d and u2d can differ on subdomain
222               igrd = 2                         ! zonal flow (m3/s) to barotropic zonal velocity (m/s)
223               DO ib = 1, SIZE(dta_alias%u2d)   ! u2d is used either over the whole bdy or only on the rim
224                  ii   = idx_bdy(jbdy)%nbi(ib,igrd)
225                  ij   = idx_bdy(jbdy)%nbj(ib,igrd)
226                  dta_alias%u2d(ib) = dta_alias%u2d(ib) / ( e2u(ii,ij) * hu_0(ii,ij) )
227               END DO
228            ENDIF
229            IF( ASSOCIATED(dta_bdy(jbdy)%v2d) ) THEN   ! no SIZE with a unassociated pointer. v2d and u2d can differ on subdomain
230               igrd = 3                         ! meridional flow (m3/s) to barotropic meridional velocity (m/s)
231               DO ib = 1, SIZE(dta_alias%v2d)   ! v2d is used either over the whole bdy or only on the rim
232                  ii   = idx_bdy(jbdy)%nbi(ib,igrd)
233                  ij   = idx_bdy(jbdy)%nbj(ib,igrd)
234                  dta_alias%v2d(ib) = dta_alias%v2d(ib) / ( e1v(ii,ij) * hv_0(ii,ij) )
235               END DO
236            ENDIF
237         ENDIF
238
239         ! tidal harmonic forcing ONLY: initialise arrays
240         IF( nn_dyn2d_dta(jbdy) == 2 ) THEN   ! we did not read ssh, u/v2d
241            IF( ASSOCIATED(dta_alias%ssh) ) dta_alias%ssh(:) = 0._wp
242            IF( ASSOCIATED(dta_alias%u2d) ) dta_alias%u2d(:) = 0._wp
243            IF( ASSOCIATED(dta_alias%v2d) ) dta_alias%v2d(:) = 0._wp
244         ENDIF
245
246         ! If full velocities in boundary data, then split it into barotropic and baroclinic component
247         IF( bf_alias(jp_bdyu3d)%ltotvel ) THEN     ! if we read 3D total velocity (can be true only if u3d was read)
248            !
249            igrd = 2                       ! zonal velocity
250            DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
251               ii   = idx_bdy(jbdy)%nbi(ib,igrd)
252               ij   = idx_bdy(jbdy)%nbj(ib,igrd)
253               dta_alias%u2d(ib) = 0._wp   ! compute barotrope zonal velocity and put it in u2d
254               DO ik = 1, jpkm1
255                  dta_alias%u2d(ib) = dta_alias%u2d(ib) + e3u_n(ii,ij,ik) * umask(ii,ij,ik) * dta_alias%u3d(ib,ik)
256               END DO
257               dta_alias%u2d(ib) =  dta_alias%u2d(ib) * r1_hu_n(ii,ij)
258               DO ik = 1, jpkm1            ! compute barocline zonal velocity and put it in u3d
259                  dta_alias%u3d(ib,ik) = dta_alias%u3d(ib,ik) - dta_alias%u2d(ib)
260               END DO
261            END DO
262            igrd = 3                       ! meridional velocity
263            DO ib = 1, idx_bdy(jbdy)%nblen(igrd)
264               ii   = idx_bdy(jbdy)%nbi(ib,igrd)
265               ij   = idx_bdy(jbdy)%nbj(ib,igrd)
266               dta_alias%v2d(ib) = 0._wp   ! compute barotrope meridional velocity and put it in v2d
267               DO ik = 1, jpkm1
268                  dta_alias%v2d(ib) = dta_alias%v2d(ib) + e3v_n(ii,ij,ik) * vmask(ii,ij,ik) * dta_alias%v3d(ib,ik)
269               END DO
270               dta_alias%v2d(ib) =  dta_alias%v2d(ib) * r1_hv_n(ii,ij)
271               DO ik = 1, jpkm1            ! compute barocline meridional velocity and put it in v3d
272                  dta_alias%v3d(ib,ik) = dta_alias%v3d(ib,ik) - dta_alias%v2d(ib)
273               END DO
274            END DO
275         ENDIF   ! ltotvel
276
277         ! update tidal harmonic forcing
278         IF( PRESENT(kit) .AND. nn_dyn2d_dta(jbdy) .GE. 2 ) THEN
279            CALL bdytide_update( kt = kt, idx = idx_bdy(jbdy), dta = dta_alias, td = tides(jbdy),   & 
280               &                 kit = kit, kt_offset = kt_offset )
281         ENDIF
282
283         !  atm surface pressure : add inverted barometer effect to ssh if it was read
284         IF ( ln_apr_obc .AND. TRIM(bf_alias(jp_bdyssh)%clrootname) /= 'NOT USED' ) THEN
285            igrd = 1
286            DO ib = 1, idx_bdy(jbdy)%nblenrim(igrd)   ! ssh is used only on the rim
287               ii   = idx_bdy(jbdy)%nbi(ib,igrd)
288               ij   = idx_bdy(jbdy)%nbj(ib,igrd)
289               dta_alias%ssh(ib) = dta_alias%ssh(ib) + ssh_ib(ii,ij)
290            END DO
291         ENDIF
292
293#if defined key_si3
294         IF( dta_alias%lneed_ice .AND. idx_bdy(jbdy)%nblen(1) > 0 ) THEN
295            ! fill temperature and salinity arrays
296            IF( TRIM(bf_alias(jp_bdyt_i)%clrootname) == 'NOT USED' )   bf_alias(jp_bdyt_i)%fnow(:,1,:) = rice_tem (jbdy)
297            IF( TRIM(bf_alias(jp_bdyt_s)%clrootname) == 'NOT USED' )   bf_alias(jp_bdyt_s)%fnow(:,1,:) = rice_tem (jbdy)
298            IF( TRIM(bf_alias(jp_bdytsu)%clrootname) == 'NOT USED' )   bf_alias(jp_bdytsu)%fnow(:,1,:) = rice_tem (jbdy)
299            IF( TRIM(bf_alias(jp_bdys_i)%clrootname) == 'NOT USED' )   bf_alias(jp_bdys_i)%fnow(:,1,:) = rice_sal (jbdy)
300            IF( TRIM(bf_alias(jp_bdyaip)%clrootname) == 'NOT USED' )   bf_alias(jp_bdyaip)%fnow(:,1,:) = rice_apnd(jbdy) * & ! rice_apnd is the pond fraction
301               &                                                                         bf_alias(jp_bdya_i)%fnow(:,1,:)     !   ( a_ip = rice_apnd * a_i )
302            IF( TRIM(bf_alias(jp_bdyhip)%clrootname) == 'NOT USED' )   bf_alias(jp_bdyhip)%fnow(:,1,:) = rice_hpnd(jbdy)
303            IF( TRIM(bf_alias(jp_bdyhil)%clrootname) == 'NOT USED' )   bf_alias(jp_bdyhil)%fnow(:,1,:) = rice_hlid(jbdy)
304
305            ! if T_i is read and not T_su, set T_su = T_i
306            IF( TRIM(bf_alias(jp_bdyt_i)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdytsu)%clrootname) == 'NOT USED' ) &
307               &   bf_alias(jp_bdytsu)%fnow(:,1,:) = bf_alias(jp_bdyt_i)%fnow(:,1,:)
308            ! if T_s is read and not T_su, set T_su = T_s
309            IF( TRIM(bf_alias(jp_bdyt_s)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdytsu)%clrootname) == 'NOT USED' ) &
310               &   bf_alias(jp_bdytsu)%fnow(:,1,:) = bf_alias(jp_bdyt_s)%fnow(:,1,:)
311            ! if T_i is read and not T_s, set T_s = T_i
312            IF( TRIM(bf_alias(jp_bdyt_i)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdyt_s)%clrootname) == 'NOT USED' ) &
313               &   bf_alias(jp_bdyt_s)%fnow(:,1,:) = bf_alias(jp_bdyt_i)%fnow(:,1,:)
314            ! if T_su is read and not T_s, set T_s = T_su
315            IF( TRIM(bf_alias(jp_bdytsu)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdyt_s)%clrootname) == 'NOT USED' ) &
316               &   bf_alias(jp_bdyt_s)%fnow(:,1,:) = bf_alias(jp_bdytsu)%fnow(:,1,:)
317            ! if T_su is read and not T_i, set T_i = (T_su + T_freeze)/2
318            IF( TRIM(bf_alias(jp_bdytsu)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdyt_i)%clrootname) == 'NOT USED' ) &
319               &   bf_alias(jp_bdyt_i)%fnow(:,1,:) = 0.5_wp * ( bf_alias(jp_bdytsu)%fnow(:,1,:) + 271.15 )
320            ! if T_s is read and not T_i, set T_i = (T_s + T_freeze)/2
321            IF( TRIM(bf_alias(jp_bdyt_s)%clrootname) /= 'NOT USED' .AND. TRIM(bf_alias(jp_bdyt_i)%clrootname) == 'NOT USED' ) &
322               &   bf_alias(jp_bdyt_i)%fnow(:,1,:) = 0.5_wp * ( bf_alias(jp_bdyt_s)%fnow(:,1,:) + 271.15 )
323
324            ! make sure ponds = 0 if no ponds scheme
325            IF ( .NOT.ln_pnd ) THEN
326               bf_alias(jp_bdyaip)%fnow(:,1,:) = 0._wp
327               bf_alias(jp_bdyhip)%fnow(:,1,:) = 0._wp
328               bf_alias(jp_bdyhil)%fnow(:,1,:) = 0._wp
329            ENDIF
330            IF ( .NOT.ln_pnd_lids ) THEN
331               bf_alias(jp_bdyhil)%fnow(:,1,:) = 0._wp
332            ENDIF
333           
334            ! convert N-cat fields (input) into jpl-cat (output)
335            ipl = SIZE(bf_alias(jp_bdya_i)%fnow, 3)           
336            IF( ipl /= jpl ) THEN      ! ice: convert N-cat fields (input) into jpl-cat (output)
337               CALL ice_var_itd( bf_alias(jp_bdyh_i)%fnow(:,1,:), bf_alias(jp_bdyh_s)%fnow(:,1,:), bf_alias(jp_bdya_i)%fnow(:,1,:), & ! in
338                  &              dta_alias%h_i                  , dta_alias%h_s                  , dta_alias%a_i                  , & ! out
339                  &              bf_alias(jp_bdyt_i)%fnow(:,1,:), bf_alias(jp_bdyt_s)%fnow(:,1,:), &                                  ! in (optional)
340                  &              bf_alias(jp_bdytsu)%fnow(:,1,:), bf_alias(jp_bdys_i)%fnow(:,1,:), &                                  ! in     -
341                  &              bf_alias(jp_bdyaip)%fnow(:,1,:), bf_alias(jp_bdyhip)%fnow(:,1,:), bf_alias(jp_bdyhil)%fnow(:,1,:), & ! in     -
342                  &              dta_alias%t_i                  , dta_alias%t_s                  , &                                  ! out    -
343                  &              dta_alias%tsu                  , dta_alias%s_i                  , &                                  ! out    -
344                  &              dta_alias%aip                  , dta_alias%hip                  , dta_alias%hil )                    ! out    -
345            ENDIF
346         ENDIF
347#endif
348      END DO  ! jbdy
349
350      IF ( ln_tide ) THEN
351         IF (ln_dynspg_ts) THEN      ! Fill temporary arrays with slow-varying bdy data                           
352            DO jbdy = 1, nb_bdy      ! Tidal component added in ts loop
353               IF ( nn_dyn2d_dta(jbdy) .GE. 2 ) THEN
354                  IF( ASSOCIATED(dta_bdy(jbdy)%ssh) ) dta_bdy_s(jbdy)%ssh(:) = dta_bdy(jbdy)%ssh(:)
355                  IF( ASSOCIATED(dta_bdy(jbdy)%u2d) ) dta_bdy_s(jbdy)%u2d(:) = dta_bdy(jbdy)%u2d(:)
356                  IF( ASSOCIATED(dta_bdy(jbdy)%v2d) ) dta_bdy_s(jbdy)%v2d(:) = dta_bdy(jbdy)%v2d(:)
357               ENDIF
358            END DO
359         ELSE ! Add tides if not split-explicit free surface else this is done in ts loop
360            !
361            CALL bdy_dta_tides( kt=kt, kt_offset=kt_offset )
362         ENDIF
363      ENDIF
364     
365      ! davbyr - add a shift to the boundary + free elevation Enda, JT from NEMO RAN 3.6
366      DO jbdy = 1, nb_bdy
367         IF( dta_bdy(jbdy)%lneed_ssh ) THEN
368            igrd  = 1
369            DO ib = 1, idx_bdy(jbdy)%nblenrim(igrd)   ! ssh is used only on the rim
370                ii = idx_bdy(jbdy)%nbi(ib,igrd)
371                ij = idx_bdy(jbdy)%nbj(ib,igrd)
372                dta_bdy(jbdy)%ssh(ib) = dta_bdy(jbdy)%ssh(ib) + rn_ssh_shift(jbdy) * tmask(ii,ij,1)
373                IF( .NOT. dta_bdy(jbdy)%lforced_ssh ) dta_bdy(jbdy)%ssh(ib) = sshn(ii,ij) * tmask(ii,ij,1)
374             END DO
375         END IF
376      END DO
377      !--- END davbyr
378     
379      !
380      IF( ln_timing )   CALL timing_stop('bdy_dta')
381      !
382   END SUBROUTINE bdy_dta
383   
384
385   SUBROUTINE bdy_dta_init
386      !!----------------------------------------------------------------------
387      !!                   ***  SUBROUTINE bdy_dta_init  ***
388      !!                   
389      !! ** Purpose :   Initialise arrays for reading of external data
390      !!                for open boundary conditions
391      !!
392      !! ** Method  :   
393      !!               
394      !!----------------------------------------------------------------------
395      INTEGER ::   jbdy, jfld    ! Local integers
396      INTEGER ::   ierror, ios     !
397      !
398      CHARACTER(len=3)                       ::   cl3           !
399      CHARACTER(len=100)                     ::   cn_dir        ! Root directory for location of data files
400      LOGICAL                                ::   ln_full_vel   ! =T => full velocities in 3D boundary data
401      !                                                         ! =F => baroclinic velocities in 3D boundary data
402      LOGICAL                                ::   ln_zinterp    ! =T => requires a vertical interpolation of the bdydta
403      REAL(wp)                               ::   rn_ice_tem, rn_ice_sal, rn_ice_age, rn_ice_apnd, rn_ice_hpnd, rn_ice_hlid
404      INTEGER                                ::   ipk,ipl       !
405      INTEGER                                ::   idvar         ! variable ID
406      INTEGER                                ::   indims        ! number of dimensions of the variable
407      INTEGER                                ::   iszdim        ! number of dimensions of the variable
408      INTEGER, DIMENSION(4)                  ::   i4dimsz       ! size of variable dimensions
409      INTEGER                                ::   igrd          ! index for grid type (1,2,3 = T,U,V)
410      LOGICAL                                ::   lluld         ! is the variable using the unlimited dimension
411      LOGICAL                                ::   llneed        !
412      LOGICAL                                ::   llread        !
413      LOGICAL                                ::   llfullbdy     !
414      TYPE(FLD_N), DIMENSION(1), TARGET  ::   bn_tem, bn_sal, bn_u3d, bn_v3d   ! must be an array to be used with fld_fill
415      TYPE(FLD_N), DIMENSION(1), TARGET  ::   bn_ssh, bn_u2d, bn_v2d           ! informations about the fields to be read
416      TYPE(FLD_N), DIMENSION(1), TARGET  ::   bn_a_i, bn_h_i, bn_h_s, bn_t_i, bn_t_s, bn_tsu, bn_s_i, bn_aip, bn_hip, bn_hil       
417      TYPE(FLD_N), DIMENSION(:), POINTER ::   bn_alias                        ! must be an array to be used with fld_fill
418      TYPE(FLD  ), DIMENSION(:), POINTER ::   bf_alias
419      !
420      NAMELIST/nambdy_dta/ cn_dir, bn_tem, bn_sal, bn_u3d, bn_v3d, bn_ssh, bn_u2d, bn_v2d,                 &
421                         & bn_a_i, bn_h_i, bn_h_s, bn_t_i, bn_t_s, bn_tsu, bn_s_i, bn_aip, bn_hip, bn_hil, &
422                         & rn_ice_tem, rn_ice_sal, rn_ice_age, rn_ice_apnd, rn_ice_hpnd, rn_ice_hlid,      &
423                         & ln_full_vel, ln_zinterp
424      !!---------------------------------------------------------------------------
425      !
426      IF(lwp) WRITE(numout,*)
427      IF(lwp) WRITE(numout,*) 'bdy_dta_ini : initialization of data at the open boundaries'
428      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
429      IF(lwp) WRITE(numout,*) ''
430
431      ALLOCATE( bf(jpbdyfld,nb_bdy), STAT=ierror )
432      IF( ierror > 0 ) THEN   
433         CALL ctl_stop( 'bdy_dta: unable to allocate bf structure' )   ;   RETURN 
434      ENDIF
435      bf(:,:)%clrootname = 'NOT USED'   ! default definition used as a flag in fld_read to do nothing.
436      bf(:,:)%lzint      = .FALSE.      ! default definition
437      bf(:,:)%ltotvel    = .FALSE.      ! default definition
438 
439      ! Read namelists
440      ! --------------
441      REWIND(numnam_cfg)
442      DO jbdy = 1, nb_bdy
443
444         WRITE(ctmp1, '(a,i2)') 'BDY number ', jbdy
445         WRITE(ctmp2, '(a,i2)') 'block nambdy_dta number ', jbdy
446
447         ! There is only one nambdy_dta block in namelist_ref -> use it for each bdy so we do a rewind
448         REWIND(numnam_ref)
449         READ  ( numnam_ref, nambdy_dta, IOSTAT = ios, ERR = 901)
450901      IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy_dta in reference namelist' )
451
452         !   by-pass nambdy_dta reading if no input data used in this bdy   
453         IF(       ( dta_bdy(jbdy)%lneed_dyn2d .AND. MOD(nn_dyn2d_dta(jbdy),2) == 1 )   &
454            & .OR. ( dta_bdy(jbdy)%lneed_dyn3d .AND.     nn_dyn3d_dta(jbdy)    == 1 )   &
455            & .OR. ( dta_bdy(jbdy)%lneed_tra   .AND.       nn_tra_dta(jbdy)    == 1 )   &
456            & .OR. ( dta_bdy(jbdy)%lneed_ice   .AND.       nn_ice_dta(jbdy)    == 1 )   )   THEN
457            ! WARNING: we don't do a rewind here, each bdy reads its own nambdy_dta block one after another
458            READ  ( numnam_cfg, nambdy_dta, IOSTAT = ios, ERR = 902 )
459902         IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy_dta in configuration namelist' )
460            IF(lwm) WRITE( numond, nambdy_dta )           
461         ENDIF
462
463         ! get the number of ice categories in bdy data file (use a_i information to do this)
464         ipl = jpl   ! default definition
465         IF( dta_bdy(jbdy)%lneed_ice ) THEN    ! if we need ice bdy data
466            IF( nn_ice_dta(jbdy) == 1 ) THEN   ! if we get ice bdy data from netcdf file
467               CALL fld_fill(  bf(jp_bdya_i,jbdy:jbdy), bn_a_i, cn_dir, 'bdy_dta', 'a_i'//' '//ctmp1, ctmp2 )   ! use namelist info
468               CALL fld_clopn( bf(jp_bdya_i,jbdy), nyear, nmonth, nday )   ! not a problem when we call it again after
469               idvar = iom_varid( bf(jp_bdya_i,jbdy)%num, bf(jp_bdya_i,jbdy)%clvar, kndims=indims, kdimsz=i4dimsz, lduld=lluld )
470               IF( indims == 4 .OR. ( indims == 3 .AND. .NOT. lluld ) ) THEN   ;   ipl = i4dimsz(3)   ! xylt or xyl
471               ELSE                                                            ;   ipl = 1            ! xy or xyt
472               ENDIF
473               bf(jp_bdya_i,jbdy)%clrootname = 'NOT USED'   ! reset to default value as this subdomain may not need to read this bdy
474            ENDIF
475         ENDIF
476
477#if defined key_si3
478         IF( .NOT.ln_pnd ) THEN
479            rn_ice_apnd = 0. ; rn_ice_hpnd = 0. ; rn_ice_hlid = 0.
480            CALL ctl_warn( 'rn_ice_apnd & rn_ice_hpnd = 0 & rn_ice_hlid = 0 when no ponds' )
481         ENDIF
482         IF( .NOT.ln_pnd_lids ) THEN
483            rn_ice_hlid = 0.
484         ENDIF
485#endif
486
487         ! temp, salt, age and ponds of incoming ice
488         rice_tem (jbdy) = rn_ice_tem
489         rice_sal (jbdy) = rn_ice_sal
490         rice_age (jbdy) = rn_ice_age
491         rice_apnd(jbdy) = rn_ice_apnd
492         rice_hpnd(jbdy) = rn_ice_hpnd
493         rice_hlid(jbdy) = rn_ice_hlid
494
495         
496         DO jfld = 1, jpbdyfld
497
498            ! =====================
499            !          ssh
500            ! =====================
501            IF( jfld == jp_bdyssh ) THEN
502               cl3 = 'ssh'
503               igrd = 1                                                    ! T point
504               ipk = 1                                                     ! surface data
505               llneed = dta_bdy(jbdy)%lneed_ssh                            ! dta_bdy(jbdy)%ssh will be needed
506               llread = MOD(nn_dyn2d_dta(jbdy),2) == 1                     ! get data from NetCDF file
507               bf_alias => bf(jp_bdyssh,jbdy:jbdy)                         ! alias for ssh structure of bdy number jbdy
508               bn_alias => bn_ssh                                          ! alias for ssh structure of nambdy_dta
509               iszdim = idx_bdy(jbdy)%nblenrim(igrd)                       ! length of this bdy on this MPI processus : used only on the rim
510            ENDIF
511            ! =====================
512            !         dyn2d
513            ! =====================
514            IF( jfld == jp_bdyu2d ) THEN
515               cl3 = 'u2d'
516               igrd = 2                                                    ! U point
517               ipk = 1                                                     ! surface data
518               llneed = dta_bdy(jbdy)%lneed_dyn2d                          ! dta_bdy(jbdy)%u2d will be needed
519               llread = .NOT. ln_full_vel .AND. MOD(nn_dyn2d_dta(jbdy),2) == 1   ! don't get u2d from u3d and read NetCDF file
520               bf_alias => bf(jp_bdyu2d,jbdy:jbdy)                         ! alias for u2d structure of bdy number jbdy
521               bn_alias => bn_u2d                                          ! alias for u2d structure of nambdy_dta
522               llfullbdy = ln_full_vel .OR. cn_dyn2d(jbdy) == 'frs'        ! need u2d over the whole bdy or only over the rim?
523               IF( llfullbdy ) THEN  ;   iszdim = idx_bdy(jbdy)%nblen(igrd)
524               ELSE                  ;   iszdim = idx_bdy(jbdy)%nblenrim(igrd)
525               ENDIF
526            ENDIF
527            IF( jfld == jp_bdyv2d ) THEN
528               cl3 = 'v2d'
529               igrd = 3                                                    ! V point
530               ipk = 1                                                     ! surface data
531               llneed = dta_bdy(jbdy)%lneed_dyn2d                          ! dta_bdy(jbdy)%v2d will be needed
532               llread = .NOT. ln_full_vel .AND. MOD(nn_dyn2d_dta(jbdy),2) == 1   ! don't get v2d from v3d and read NetCDF file
533               bf_alias => bf(jp_bdyv2d,jbdy:jbdy)                         ! alias for v2d structure of bdy number jbdy
534               bn_alias => bn_v2d                                          ! alias for v2d structure of nambdy_dta
535               llfullbdy = ln_full_vel .OR. cn_dyn2d(jbdy) == 'frs'        ! need v2d over the whole bdy or only over the rim?
536               IF( llfullbdy ) THEN  ;   iszdim = idx_bdy(jbdy)%nblen(igrd)
537               ELSE                  ;   iszdim = idx_bdy(jbdy)%nblenrim(igrd)
538               ENDIF
539            ENDIF
540            ! =====================
541            !         dyn3d
542            ! =====================
543            IF( jfld == jp_bdyu3d ) THEN
544               cl3 = 'u3d'
545               igrd = 2                                                    ! U point
546               ipk = jpk                                                   ! 3d data
547               llneed = dta_bdy(jbdy)%lneed_dyn3d .OR.               &     ! dta_bdy(jbdy)%u3d will be needed
548                  &   ( dta_bdy(jbdy)%lneed_dyn2d .AND. ln_full_vel )      !   u3d needed to compute u2d
549               llread = nn_dyn3d_dta(jbdy) == 1                            ! get data from NetCDF file
550               bf_alias => bf(jp_bdyu3d,jbdy:jbdy)                         ! alias for u3d structure of bdy number jbdy
551               bn_alias => bn_u3d                                          ! alias for u3d structure of nambdy_dta
552               iszdim = idx_bdy(jbdy)%nblen(igrd)                          ! length of this bdy on this MPI processus
553           ENDIF
554            IF( jfld == jp_bdyv3d ) THEN
555               cl3 = 'v3d'
556               igrd = 3                                                    ! V point
557               ipk = jpk                                                   ! 3d data
558               llneed = dta_bdy(jbdy)%lneed_dyn3d .OR.               &     ! dta_bdy(jbdy)%v3d will be needed
559                  &   ( dta_bdy(jbdy)%lneed_dyn2d .AND. ln_full_vel )      !   v3d needed to compute v2d
560               llread = nn_dyn3d_dta(jbdy) == 1                            ! get data from NetCDF file
561               bf_alias => bf(jp_bdyv3d,jbdy:jbdy)                         ! alias for v3d structure of bdy number jbdy
562               bn_alias => bn_v3d                                          ! alias for v3d structure of nambdy_dta
563               iszdim = idx_bdy(jbdy)%nblen(igrd)                          ! length of this bdy on this MPI processus
564           ENDIF
565
566            ! =====================
567            !          tra
568            ! =====================
569            IF( jfld == jp_bdytem ) THEN
570               cl3 = 'tem'
571               igrd = 1                                                    ! T point
572               ipk = jpk                                                   ! 3d data
573               llneed = dta_bdy(jbdy)%lneed_tra                            ! dta_bdy(jbdy)%tem will be needed
574               llread = nn_tra_dta(jbdy) == 1                              ! get data from NetCDF file
575               bf_alias => bf(jp_bdytem,jbdy:jbdy)                         ! alias for ssh structure of bdy number jbdy
576               bn_alias => bn_tem                                          ! alias for ssh structure of nambdy_dta
577               iszdim = idx_bdy(jbdy)%nblen(igrd)                          ! length of this bdy on this MPI processus
578            ENDIF
579            IF( jfld == jp_bdysal ) THEN
580               cl3 = 'sal'
581               igrd = 1                                                    ! T point
582               ipk = jpk                                                   ! 3d data
583               llneed = dta_bdy(jbdy)%lneed_tra                            ! dta_bdy(jbdy)%sal will be needed
584               llread = nn_tra_dta(jbdy) == 1                              ! get data from NetCDF file
585               bf_alias => bf(jp_bdysal,jbdy:jbdy)                         ! alias for ssh structure of bdy number jbdy
586               bn_alias => bn_sal                                          ! alias for ssh structure of nambdy_dta
587               iszdim = idx_bdy(jbdy)%nblen(igrd)                          ! length of this bdy on this MPI processus
588            ENDIF
589
590            ! =====================
591            !          ice
592            ! =====================
593            IF(  jfld == jp_bdya_i .OR. jfld == jp_bdyh_i .OR. jfld == jp_bdyh_s .OR. &
594               & jfld == jp_bdyt_i .OR. jfld == jp_bdyt_s .OR. jfld == jp_bdytsu .OR. &
595               & jfld == jp_bdys_i .OR. jfld == jp_bdyaip .OR. jfld == jp_bdyhip .OR. jfld == jp_bdyhil ) THEN
596               igrd = 1                                                    ! T point
597               ipk = ipl                                                   ! jpl-cat data
598               llneed = dta_bdy(jbdy)%lneed_ice                            ! ice will be needed
599               llread = nn_ice_dta(jbdy) == 1                              ! get data from NetCDF file
600               iszdim = idx_bdy(jbdy)%nblen(igrd)                          ! length of this bdy on this MPI processus
601            ENDIF
602            IF( jfld == jp_bdya_i ) THEN
603               cl3 = 'a_i'
604               bf_alias => bf(jp_bdya_i,jbdy:jbdy)                         ! alias for a_i structure of bdy number jbdy
605               bn_alias => bn_a_i                                          ! alias for a_i structure of nambdy_dta
606            ENDIF
607            IF( jfld == jp_bdyh_i ) THEN
608               cl3 = 'h_i'
609               bf_alias => bf(jp_bdyh_i,jbdy:jbdy)                         ! alias for h_i structure of bdy number jbdy
610               bn_alias => bn_h_i                                          ! alias for h_i structure of nambdy_dta
611            ENDIF
612            IF( jfld == jp_bdyh_s ) THEN
613               cl3 = 'h_s'
614               bf_alias => bf(jp_bdyh_s,jbdy:jbdy)                         ! alias for h_s structure of bdy number jbdy
615               bn_alias => bn_h_s                                          ! alias for h_s structure of nambdy_dta
616            ENDIF
617            IF( jfld == jp_bdyt_i ) THEN
618               cl3 = 't_i'
619               bf_alias => bf(jp_bdyt_i,jbdy:jbdy)                         ! alias for t_i structure of bdy number jbdy
620               bn_alias => bn_t_i                                          ! alias for t_i structure of nambdy_dta
621            ENDIF
622            IF( jfld == jp_bdyt_s ) THEN
623               cl3 = 't_s'
624               bf_alias => bf(jp_bdyt_s,jbdy:jbdy)                         ! alias for t_s structure of bdy number jbdy
625               bn_alias => bn_t_s                                          ! alias for t_s structure of nambdy_dta
626            ENDIF
627            IF( jfld == jp_bdytsu ) THEN
628               cl3 = 'tsu'
629               bf_alias => bf(jp_bdytsu,jbdy:jbdy)                         ! alias for tsu structure of bdy number jbdy
630               bn_alias => bn_tsu                                          ! alias for tsu structure of nambdy_dta
631            ENDIF
632            IF( jfld == jp_bdys_i ) THEN
633               cl3 = 's_i'
634               bf_alias => bf(jp_bdys_i,jbdy:jbdy)                         ! alias for s_i structure of bdy number jbdy
635               bn_alias => bn_s_i                                          ! alias for s_i structure of nambdy_dta
636            ENDIF
637            IF( jfld == jp_bdyaip ) THEN
638               cl3 = 'aip'
639               bf_alias => bf(jp_bdyaip,jbdy:jbdy)                         ! alias for aip structure of bdy number jbdy
640               bn_alias => bn_aip                                          ! alias for aip structure of nambdy_dta
641            ENDIF
642            IF( jfld == jp_bdyhip ) THEN
643               cl3 = 'hip'
644               bf_alias => bf(jp_bdyhip,jbdy:jbdy)                         ! alias for hip structure of bdy number jbdy
645               bn_alias => bn_hip                                          ! alias for hip structure of nambdy_dta
646            ENDIF
647            IF( jfld == jp_bdyhil ) THEN
648               cl3 = 'hil'
649               bf_alias => bf(jp_bdyhil,jbdy:jbdy)                         ! alias for hil structure of bdy number jbdy
650               bn_alias => bn_hil                                          ! alias for hil structure of nambdy_dta
651            ENDIF
652
653            IF( llneed .AND. iszdim > 0 ) THEN                             ! dta_bdy(jbdy)%xxx will be needed
654               !                                                           !   -> must be associated with an allocated target
655               ALLOCATE( bf_alias(1)%fnow( iszdim, 1, ipk ) )              ! allocate the target
656               !
657               IF( llread ) THEN                                           ! get data from NetCDF file
658                  CALL fld_fill( bf_alias, bn_alias, cn_dir, 'bdy_dta', cl3//' '//ctmp1, ctmp2 )   ! use namelist info
659                  IF( bf_alias(1)%ln_tint ) ALLOCATE( bf_alias(1)%fdta( iszdim, 1, ipk, 2 ) )
660                  bf_alias(1)%imap    => idx_bdy(jbdy)%nbmap(1:iszdim,igrd)   ! associate the mapping used for this bdy
661                  bf_alias(1)%igrd    = igrd                                  ! used only for vertical integration of 3D arrays
662                  bf_alias(1)%ibdy    = jbdy                                  !  "    "    "     "          "      "  "    "   
663                  bf_alias(1)%ltotvel = ln_full_vel                           ! T if u3d is full velocity
664                  bf_alias(1)%lzint   = ln_zinterp                            ! T if it requires a vertical interpolation
665               ENDIF
666
667               ! associate the pointer and get rid of the dimensions with a size equal to 1
668               IF( jfld == jp_bdyssh )        dta_bdy(jbdy)%ssh => bf_alias(1)%fnow(:,1,1)
669               IF( jfld == jp_bdyu2d )        dta_bdy(jbdy)%u2d => bf_alias(1)%fnow(:,1,1)
670               IF( jfld == jp_bdyv2d )        dta_bdy(jbdy)%v2d => bf_alias(1)%fnow(:,1,1)
671               IF( jfld == jp_bdyu3d )        dta_bdy(jbdy)%u3d => bf_alias(1)%fnow(:,1,:)
672               IF( jfld == jp_bdyv3d )        dta_bdy(jbdy)%v3d => bf_alias(1)%fnow(:,1,:)
673               IF( jfld == jp_bdytem )        dta_bdy(jbdy)%tem => bf_alias(1)%fnow(:,1,:)
674               IF( jfld == jp_bdysal )        dta_bdy(jbdy)%sal => bf_alias(1)%fnow(:,1,:)
675               IF( jfld == jp_bdya_i ) THEN
676                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%a_i => bf_alias(1)%fnow(:,1,:)
677                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%a_i(iszdim,jpl) )
678                  ENDIF
679               ENDIF
680               IF( jfld == jp_bdyh_i ) THEN
681                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%h_i => bf_alias(1)%fnow(:,1,:)
682                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%h_i(iszdim,jpl) )
683                  ENDIF
684               ENDIF
685               IF( jfld == jp_bdyh_s ) THEN
686                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%h_s => bf_alias(1)%fnow(:,1,:)
687                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%h_s(iszdim,jpl) )
688                  ENDIF
689               ENDIF
690               IF( jfld == jp_bdyt_i ) THEN
691                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%t_i => bf_alias(1)%fnow(:,1,:)
692                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%t_i(iszdim,jpl) )
693                  ENDIF
694               ENDIF
695               IF( jfld == jp_bdyt_s ) THEN
696                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%t_s => bf_alias(1)%fnow(:,1,:)
697                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%t_s(iszdim,jpl) )
698                  ENDIF
699               ENDIF
700               IF( jfld == jp_bdytsu ) THEN
701                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%tsu => bf_alias(1)%fnow(:,1,:)
702                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%tsu(iszdim,jpl) )
703                  ENDIF
704               ENDIF
705               IF( jfld == jp_bdys_i ) THEN
706                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%s_i => bf_alias(1)%fnow(:,1,:)
707                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%s_i(iszdim,jpl) )
708                  ENDIF
709               ENDIF
710               IF( jfld == jp_bdyaip ) THEN
711                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%aip => bf_alias(1)%fnow(:,1,:)
712                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%aip(iszdim,jpl) )
713                  ENDIF
714               ENDIF
715               IF( jfld == jp_bdyhip ) THEN
716                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%hip => bf_alias(1)%fnow(:,1,:)
717                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%hip(iszdim,jpl) )
718                  ENDIF
719               ENDIF
720               IF( jfld == jp_bdyhil ) THEN
721                  IF( ipk == jpl ) THEN   ;   dta_bdy(jbdy)%hil => bf_alias(1)%fnow(:,1,:)
722                  ELSE                    ;   ALLOCATE( dta_bdy(jbdy)%hil(iszdim,jpl) )
723                  ENDIF
724               ENDIF
725            ENDIF
726
727         END DO   ! jpbdyfld
728         !
729      END DO ! jbdy
730      !
731   END SUBROUTINE bdy_dta_init
732   
733   !!==============================================================================
734END MODULE bdydta
Note: See TracBrowser for help on using the repository browser.