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 branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/BDY/bdydta.F90 @ 4292

Last change on this file since 4292 was 4292, checked in by cetlod, 10 years ago

dev_MERGE_2013 : 1st step of the merge, see ticket #1185

  • Property svn:keywords set to Id
File size: 40.2 KB
RevLine 
[911]1MODULE bdydta
[1125]2   !!======================================================================
3   !!                       ***  MODULE bdydta  ***
[911]4   !! Open boundary data : read the data for the unstructured open boundaries.
[1125]5   !!======================================================================
6   !! History :  1.0  !  2005-01  (J. Chanut, A. Sellar)  Original code
7   !!             -   !  2007-01  (D. Storkey) Update to use IOM module
[2528]8   !!             -   !  2007-07  (D. Storkey) add bdy_dta_fla
[1125]9   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
[2528]10   !!            3.3  !  2010-09  (E.O'Dea) modifications for Shelf configurations
11   !!            3.3  !  2010-09  (D.Storkey) add ice boundary conditions
[3294]12   !!            3.4  !  2011     (D. Storkey) rewrite in preparation for OBC-BDY merge
[4292]13   !!            3.6  !  2012-01  (C. Rousset) add ice boundary conditions for lim3
[1125]14   !!----------------------------------------------------------------------
15#if defined key_bdy
16   !!----------------------------------------------------------------------
[3294]17   !!   'key_bdy'                     Open Boundary Conditions
[1125]18   !!----------------------------------------------------------------------
[3294]19   !!    bdy_dta        : read external data along open boundaries from file
20   !!    bdy_dta_init   : initialise arrays etc for reading of external data
[1125]21   !!----------------------------------------------------------------------
[3294]22   USE wrk_nemo        ! Memory Allocation
23   USE timing          ! Timing
[911]24   USE oce             ! ocean dynamics and tracers
25   USE dom_oce         ! ocean space and time domain
26   USE phycst          ! physical constants
[3294]27   USE bdy_oce         ! ocean open boundary conditions 
[911]28   USE bdytides        ! tidal forcing at boundaries
[3294]29   USE fldread         ! read input fields
30   USE iom             ! IOM library
[911]31   USE in_out_manager  ! I/O logical units
[4292]32   USE dynspg_oce, ONLY: lk_dynspg_ts ! Split-explicit free surface flag
[2528]33#if defined key_lim2
34   USE ice_2
[4292]35#elif defined key_lim3
36   USE par_ice
37   USE ice
38   USE limcat_1D          ! redistribute ice input into categories
[2528]39#endif
[3651]40   USE sbcapr
[911]41
42   IMPLICIT NONE
43   PRIVATE
44
[3294]45   PUBLIC   bdy_dta          ! routine called by step.F90 and dynspg_ts.F90
46   PUBLIC   bdy_dta_init     ! routine called by nemogcm.F90
[911]47
[3294]48   INTEGER, ALLOCATABLE, DIMENSION(:)   ::   nb_bdy_fld        ! Number of fields to update for each boundary set.
49   INTEGER                              ::   nb_bdy_fld_sum    ! Total number of fields to update for all boundary sets.
[911]50
[3294]51   LOGICAL,           DIMENSION(jp_bdy) ::   ln_full_vel_array ! =T => full velocities in 3D boundary conditions
52                                                               ! =F => baroclinic velocities in 3D boundary conditions
[911]53
[3294]54   TYPE(FLD), PUBLIC, ALLOCATABLE, DIMENSION(:), TARGET ::   bf        ! structure of input fields (file informations, fields read)
[911]55
[3294]56   TYPE(MAP_POINTER), ALLOCATABLE, DIMENSION(:) :: nbmap_ptr   ! array of pointers to nbmap
[911]57
[4292]58#if defined key_lim3
59   LOGICAL :: ll_bdylim3                  ! determine whether ice input is lim2 (F) or lim3 (T) type
60   INTEGER :: jfld_hti, jfld_hts, jfld_ai ! indices of ice thickness, snow thickness and concentration in bf structure
61#endif
62
[3294]63#  include "domzgr_substitute.h90"
[1125]64   !!----------------------------------------------------------------------
[2528]65   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[1146]66   !! $Id$
[2528]67   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
[1125]68   !!----------------------------------------------------------------------
[911]69CONTAINS
70
[3294]71      SUBROUTINE bdy_dta( kt, jit, time_offset )
[1125]72      !!----------------------------------------------------------------------
[3294]73      !!                   ***  SUBROUTINE bdy_dta  ***
[911]74      !!                   
[3294]75      !! ** Purpose :   Update external data for open boundary conditions
[911]76      !!
[3294]77      !! ** Method  :   Use fldread.F90
78      !!               
[1125]79      !!----------------------------------------------------------------------
[911]80      !!
[3294]81      INTEGER, INTENT( in )           ::   kt    ! ocean time-step index
82      INTEGER, INTENT( in ), OPTIONAL ::   jit   ! subcycle time-step index (for timesplitting option)
83      INTEGER, INTENT( in ), OPTIONAL ::   time_offset  ! time offset in units of timesteps. NB. if jit
84                                                        ! is present then units = subcycle timesteps.
85                                                        ! time_offset = 0 => get data at "now" time level
86                                                        ! time_offset = -1 => get data at "before" time level
87                                                        ! time_offset = +1 => get data at "after" time level
88                                                        ! etc.
89      !!
[4292]90      INTEGER     ::  ib_bdy, jfld, jstart, jend, ib, ii, ij, ik, igrd, jl  ! local indices
[3294]91      INTEGER,          DIMENSION(jpbgrd) ::   ilen1 
92      INTEGER, POINTER, DIMENSION(:)      ::   nblen, nblenrim  ! short cuts
[4292]93      TYPE(OBC_DATA), POINTER             ::   dta              ! short cut
[3294]94      !!
[911]95      !!---------------------------------------------------------------------------
[3294]96      !!
97      IF( nn_timing == 1 ) CALL timing_start('bdy_dta')
[911]98
[3294]99      ! Initialise data arrays once for all from initial conditions where required
100      !---------------------------------------------------------------------------
101      IF( kt .eq. nit000 .and. .not. PRESENT(jit) ) THEN
[1125]102
[3294]103         ! Calculate depth-mean currents
104         !-----------------------------
[4292]105         CALL wrk_alloc(jpi,jpj,pun2d,pvn2d) 
[2528]106
[4292]107         pun2d(:,:) = 0.e0
108         pvn2d(:,:) = 0.e0
[3294]109         DO ik = 1, jpkm1   !! Vertically integrated momentum trends
[4292]110             pun2d(:,:) = pun2d(:,:) + fse3u(:,:,ik) * umask(:,:,ik) * un(:,:,ik)
111             pvn2d(:,:) = pvn2d(:,:) + fse3v(:,:,ik) * vmask(:,:,ik) * vn(:,:,ik)
[3294]112         END DO
[4292]113         pun2d(:,:) = pun2d(:,:) * hur(:,:)
114         pvn2d(:,:) = pvn2d(:,:) * hvr(:,:)
[3294]115         
116         DO ib_bdy = 1, nb_bdy
[911]117
[3294]118            nblen => idx_bdy(ib_bdy)%nblen
119            nblenrim => idx_bdy(ib_bdy)%nblenrim
[4292]120            dta => dta_bdy(ib_bdy)
[2528]121
[4292]122            IF( nn_dyn2d_dta(ib_bdy) .eq. 0 ) THEN
[3651]123               ilen1(:) = nblen(:)
[4292]124               IF( dta%ll_ssh ) THEN
125                  igrd = 1
126                  DO ib = 1, ilen1(igrd)
[3294]127                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
128                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
[4292]129                     dta_bdy(ib_bdy)%ssh(ib) = sshn(ii,ij) * tmask(ii,ij,1)         
130                  END DO
131               END IF
132               IF( dta%ll_u2d ) THEN
133                  igrd = 2
134                  DO ib = 1, ilen1(igrd)
[3294]135                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
136                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
[4292]137                     dta_bdy(ib_bdy)%u2d(ib) = pun2d(ii,ij) * umask(ii,ij,1)         
138                  END DO
139               END IF
140               IF( dta%ll_v2d ) THEN
141                  igrd = 3
142                  DO ib = 1, ilen1(igrd)
143                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
144                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
145                     dta_bdy(ib_bdy)%v2d(ib) = pvn2d(ii,ij) * vmask(ii,ij,1)         
146                  END DO
147               END IF
148            ENDIF
149
150            IF( nn_dyn3d_dta(ib_bdy) .eq. 0 ) THEN
151               ilen1(:) = nblen(:)
152               IF( dta%ll_u3d ) THEN
153                  igrd = 2 
154                  DO ib = 1, ilen1(igrd)
155                     DO ik = 1, jpkm1
156                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
157                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
158                        dta_bdy(ib_bdy)%u3d(ib,ik) =  ( un(ii,ij,ik) - pun2d(ii,ij) ) * umask(ii,ij,ik)         
[3294]159                     END DO
[4292]160                  END DO
161               END IF
162               IF( dta%ll_v3d ) THEN
163                  igrd = 3 
164                  DO ib = 1, ilen1(igrd)
165                     DO ik = 1, jpkm1
166                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
167                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
168                        dta_bdy(ib_bdy)%v3d(ib,ik) =  ( vn(ii,ij,ik) - pvn2d(ii,ij) ) * vmask(ii,ij,ik)         
169                        END DO
170                  END DO
171               END IF
[3294]172            ENDIF
[911]173
[4292]174            IF( nn_tra_dta(ib_bdy) .eq. 0 ) THEN
[3651]175               ilen1(:) = nblen(:)
[4292]176               IF( dta%ll_tem ) THEN
177                  igrd = 1 
178                  DO ib = 1, ilen1(igrd)
179                     DO ik = 1, jpkm1
180                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
181                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
182                        dta_bdy(ib_bdy)%tem(ib,ik) = tsn(ii,ij,ik,jp_tem) * tmask(ii,ij,ik)         
183                     END DO
184                  END DO
185               END IF
186               IF( dta%ll_sal ) THEN
187                  igrd = 1 
188                  DO ib = 1, ilen1(igrd)
189                     DO ik = 1, jpkm1
190                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
191                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
192                        dta_bdy(ib_bdy)%sal(ib,ik) = tsn(ii,ij,ik,jp_sal) * tmask(ii,ij,ik)         
193                     END DO
194                  END DO
195               END IF
[3294]196            ENDIF
[911]197
[3294]198#if defined key_lim2
[4292]199            IF( nn_ice_lim2_dta(ib_bdy) .eq. 0 ) THEN
[3651]200               ilen1(:) = nblen(:)
[4292]201               IF( dta%ll_frld ) THEN
202                  igrd = 1 
203                  DO ib = 1, ilen1(igrd)
204                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
205                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
206                     dta_bdy(ib_bdy)%frld(ib) = frld(ii,ij) * tmask(ii,ij,1)         
207                  END DO
208               END IF
209               IF( dta%ll_hicif ) THEN
210                  igrd = 1 
211                  DO ib = 1, ilen1(igrd)
212                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
213                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
214                     dta_bdy(ib_bdy)%hicif(ib) = hicif(ii,ij) * tmask(ii,ij,1)         
215                  END DO
216               END IF
217               IF( dta%ll_hsnif ) THEN
218                  igrd = 1 
219                  DO ib = 1, ilen1(igrd)
220                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
221                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
222                     dta_bdy(ib_bdy)%hsnif(ib) = hsnif(ii,ij) * tmask(ii,ij,1)         
223                  END DO
224               END IF
[1125]225            ENDIF
[4292]226#elif defined key_lim3
227            IF( nn_ice_lim_dta(ib_bdy) .eq. 0 ) THEN
228               ilen1(:) = nblen(:)
229               IF( dta%ll_a_i ) THEN
230                  igrd = 1   
231                  DO jl = 1, jpl
232                     DO ib = 1, ilen1(igrd)
233                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
234                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
235                        dta_bdy(ib_bdy)%a_i (ib,jl) =  a_i(ii,ij,jl) * tmask(ii,ij,1) 
236                     END DO
237                  END DO
238               ENDIF
239               IF( dta%ll_ht_i ) THEN
240                  igrd = 1   
241                  DO jl = 1, jpl
242                     DO ib = 1, ilen1(igrd)
243                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
244                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
245                        dta_bdy(ib_bdy)%ht_i (ib,jl) =  ht_i(ii,ij,jl) * tmask(ii,ij,1) 
246                     END DO
247                  END DO
248               ENDIF
249               IF( dta%ll_ht_s ) THEN
250                  igrd = 1   
251                  DO jl = 1, jpl
252                     DO ib = 1, ilen1(igrd)
253                        ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
254                        ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
255                        dta_bdy(ib_bdy)%ht_s (ib,jl) =  ht_s(ii,ij,jl) * tmask(ii,ij,1) 
256                     END DO
257                  END DO
258               ENDIF
259            ENDIF
[3294]260#endif
[911]261
[3294]262         ENDDO ! ib_bdy
[911]263
[4292]264         CALL wrk_dealloc(jpi,jpj,pun2d,pvn2d) 
[911]265
[3294]266      ENDIF ! kt .eq. nit000
[911]267
[3294]268      ! update external data from files
269      !--------------------------------
270     
271      jstart = 1
272      DO ib_bdy = 1, nb_bdy   
[4292]273         dta => dta_bdy(ib_bdy)
[3294]274         IF( nn_dta(ib_bdy) .eq. 1 ) THEN ! skip this bit if no external data required
275     
276            IF( PRESENT(jit) ) THEN
277               ! Update barotropic boundary conditions only
[3651]278               ! jit is optional argument for fld_read and bdytide_update
[4292]279               IF( cn_dyn2d(ib_bdy) /= 'none' ) THEN
[3294]280                  IF( nn_dyn2d_dta(ib_bdy) .eq. 2 ) THEN ! tidal harmonic forcing ONLY: initialise arrays
[4292]281                     IF( dta%ll_ssh ) dta%ssh(:) = 0.0
282                     IF( dta%ll_u2d ) dta%u2d(:) = 0.0
283                     IF( dta%ll_u3d ) dta%v2d(:) = 0.0
[3294]284                  ENDIF
[4292]285                  IF (cn_tra(ib_bdy) /= 'runoff') THEN
286                     IF( nn_dyn2d_dta(ib_bdy) .EQ. 1 .OR. nn_dyn2d_dta(ib_bdy) .EQ. 3 ) THEN
[3703]287
[4292]288                        jend = jstart + dta%nread(2) - 1
[3703]289                        CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), map=nbmap_ptr(jstart:jend),  &
[3851]290                                     & kit=jit, kt_offset=time_offset )
[3703]291
[4292]292                        ! If full velocities in boundary data then extract barotropic velocities from 3D fields
[3703]293                        IF( ln_full_vel_array(ib_bdy) .AND.                                             &
294                          &    ( nn_dyn2d_dta(ib_bdy) .EQ. 1 .OR. nn_dyn2d_dta(ib_bdy) .EQ. 3 .OR.  &
295                          &      nn_dyn3d_dta(ib_bdy) .EQ. 1 ) )THEN
296
[3651]297                           igrd = 2                      ! zonal velocity
[4292]298                           dta%u2d(:) = 0.0
[3651]299                           DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
300                              ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
301                              ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
302                              DO ik = 1, jpkm1
[4292]303                                 dta%u2d(ib) = dta%u2d(ib) &
304                       &                          + fse3u(ii,ij,ik) * umask(ii,ij,ik) * dta%u3d(ib,ik)
[3651]305                              END DO
[4292]306                              dta%u2d(ib) =  dta%u2d(ib) * hur(ii,ij)
[3651]307                           END DO
308                           igrd = 3                      ! meridional velocity
[4292]309                           dta%v2d(:) = 0.0
[3651]310                           DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
311                              ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
312                              ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
313                              DO ik = 1, jpkm1
[4292]314                                 dta%v2d(ib) = dta%v2d(ib) &
315                       &                       + fse3v(ii,ij,ik) * vmask(ii,ij,ik) * dta%v3d(ib,ik)
[3651]316                              END DO
[4292]317                              dta%v2d(ib) =  dta%v2d(ib) * hvr(ii,ij)
[3651]318                           END DO
319                        ENDIF                   
320                     ENDIF
321                     IF( nn_dyn2d_dta(ib_bdy) .ge. 2 ) THEN ! update tidal harmonic forcing
[4292]322                        CALL bdytide_update( kt=kt, idx=idx_bdy(ib_bdy), dta=dta, td=tides(ib_bdy),   & 
[3651]323                          &                 jit=jit, time_offset=time_offset )
324                     ENDIF
[3294]325                  ENDIF
[1125]326               ENDIF
[3294]327            ELSE
[4292]328               IF (cn_tra(ib_bdy) == 'runoff') then      ! runoff condition
[3651]329                  jend = nb_bdy_fld(ib_bdy)
[3703]330                  CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend),  &
[3851]331                               & map=nbmap_ptr(jstart:jend), kt_offset=time_offset )
[3651]332                  !
333                  igrd = 2                      ! zonal velocity
334                  DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
335                     ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
336                     ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
[4292]337                     dta%u2d(ib) = dta%u2d(ib) / ( e2u(ii,ij) * hu_0(ii,ij) )
[3651]338                  END DO
339                  !
340                  igrd = 3                      ! meridional velocity
341                  DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
342                     ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
343                     ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
[4292]344                     dta%v2d(ib) = dta%v2d(ib) / ( e1v(ii,ij) * hv_0(ii,ij) )
[3651]345                  END DO
346               ELSE
[4292]347                  IF( nn_dyn2d_dta(ib_bdy) .eq. 2 ) THEN ! tidal harmonic forcing ONLY: initialise arrays
348                     IF( dta%ll_ssh ) dta%ssh(:) = 0.0
349                     IF( dta%ll_u2d ) dta%u2d(:) = 0.0
350                     IF( dta%ll_v2d ) dta%v2d(:) = 0.0
[3651]351                  ENDIF
[4292]352                  IF( dta%nread(1) .gt. 0 ) THEN ! update external data
353                     jend = jstart + dta%nread(1) - 1
[3703]354                     CALL fld_read( kt=kt, kn_fsbc=1, sd=bf(jstart:jend), &
[3851]355                                  & map=nbmap_ptr(jstart:jend), kt_offset=time_offset )
[3651]356                  ENDIF
357                  ! If full velocities in boundary data then split into barotropic and baroclinic data
358                  IF( ln_full_vel_array(ib_bdy) .and.                                             &
[3703]359                    & ( nn_dyn2d_dta(ib_bdy) .EQ. 1 .OR. nn_dyn2d_dta(ib_bdy) .EQ. 3 .OR. &
360                    &   nn_dyn3d_dta(ib_bdy) .EQ. 1 ) ) THEN
[3651]361                     igrd = 2                      ! zonal velocity
[4292]362                     dta%u2d(:) = 0.0
[3651]363                     DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
364                        ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
365                        ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
366                        DO ik = 1, jpkm1
[4292]367                           dta%u2d(ib) = dta%u2d(ib) &
368                 &                       + fse3u(ii,ij,ik) * umask(ii,ij,ik) * dta%u3d(ib,ik)
[3651]369                        END DO
[4292]370                        dta%u2d(ib) =  dta%u2d(ib) * hur(ii,ij)
[3651]371                        DO ik = 1, jpkm1
[4292]372                           dta%u3d(ib,ik) = dta%u3d(ib,ik) - dta%u2d(ib)
[3651]373                        END DO
374                     END DO
375                     igrd = 3                      ! meridional velocity
[4292]376                     dta%v2d(:) = 0.0
[3651]377                     DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
378                        ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
379                        ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
380                        DO ik = 1, jpkm1
[4292]381                           dta%v2d(ib) = dta%v2d(ib) &
382                 &                       + fse3v(ii,ij,ik) * vmask(ii,ij,ik) * dta%v3d(ib,ik)
[3651]383                        END DO
[4292]384                        dta%v2d(ib) =  dta%v2d(ib) * hvr(ii,ij)
[3651]385                        DO ik = 1, jpkm1
[4292]386                           dta%v3d(ib,ik) = dta%v3d(ib,ik) - dta%v2d(ib)
[3651]387                        END DO
388                     END DO
389                  ENDIF
[4292]390
[1125]391               ENDIF
[4292]392#if defined key_lim3
393               IF( .NOT. ll_bdylim3 .AND. nn_ice_lim(ib_bdy) > 0 .AND. nn_ice_lim_dta(ib_bdy) == 1 ) THEN ! bdy ice input (case input is lim2 type)
394                CALL lim_cat_1D ( bf(jfld_hti)%fnow(:,1,1), bf(jfld_hts)%fnow(:,1,1), bf(jfld_ai)%fnow(:,1,1), &
395                                  & dta_bdy(ib_bdy)%ht_i,     dta_bdy(ib_bdy)%ht_s,     dta_bdy(ib_bdy)%a_i     )
396               ENDIF
397#endif
[1125]398            ENDIF
[4292]399            jstart = jstart + dta%nread(1)
[3651]400         END IF ! nn_dta(ib_bdy) = 1
401      END DO  ! ib_bdy
[911]402
[4292]403      ! bg jchanut tschanges
404#if defined key_tide
405      ! Add tides if not split-explicit free surface else this is done in ts loop
406      IF (.NOT.lk_dynspg_ts) CALL bdy_dta_tides( kt=kt, time_offset=time_offset )
407#endif
408      ! end jchanut tschanges
409
[3651]410      IF ( ln_apr_obc ) THEN
411         DO ib_bdy = 1, nb_bdy
[4292]412            IF (cn_tra(ib_bdy) /= 'runoff')THEN
[3651]413               igrd = 1                      ! meridional velocity
414               DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
[3294]415                  ii   = idx_bdy(ib_bdy)%nbi(ib,igrd)
416                  ij   = idx_bdy(ib_bdy)%nbj(ib,igrd)
[3651]417                  dta_bdy(ib_bdy)%ssh(ib) = dta_bdy(ib_bdy)%ssh(ib) + ssh_ib(ii,ij)
418               ENDDO
[1125]419            ENDIF
[3651]420         ENDDO
421      ENDIF
[911]422
[3294]423      IF( nn_timing == 1 ) CALL timing_stop('bdy_dta')
[911]424
[3294]425      END SUBROUTINE bdy_dta
[911]426
427
[3294]428      SUBROUTINE bdy_dta_init
429      !!----------------------------------------------------------------------
430      !!                   ***  SUBROUTINE bdy_dta_init  ***
431      !!                   
432      !! ** Purpose :   Initialise arrays for reading of external data
433      !!                for open boundary conditions
434      !!
[4292]435      !! ** Method  :   
[3294]436      !!               
437      !!----------------------------------------------------------------------
438      USE dynspg_oce, ONLY: lk_dynspg_ts
439      !!
440      INTEGER     ::  ib_bdy, jfld, jstart, jend, ierror  ! local indices
[4147]441      INTEGER      ::   ios                               ! Local integer output status for namelist read
[3294]442      !!
443      CHARACTER(len=100)                     ::   cn_dir        ! Root directory for location of data files
444      CHARACTER(len=100), DIMENSION(nb_bdy)  ::   cn_dir_array  ! Root directory for location of data files
445      LOGICAL                                ::   ln_full_vel   ! =T => full velocities in 3D boundary data
446                                                                ! =F => baroclinic velocities in 3D boundary data
447      INTEGER                                ::   ilen_global   ! Max length required for global bdy dta arrays
448      INTEGER, ALLOCATABLE, DIMENSION(:)     ::   ilen1, ilen3  ! size of 1st and 3rd dimensions of local arrays
449      INTEGER, ALLOCATABLE, DIMENSION(:)     ::   ibdy           ! bdy set for a particular jfld
450      INTEGER, ALLOCATABLE, DIMENSION(:)     ::   igrid         ! index for grid type (1,2,3 = T,U,V)
451      INTEGER, POINTER, DIMENSION(:)         ::   nblen, nblenrim  ! short cuts
[4292]452      TYPE(OBC_DATA), POINTER                ::   dta           ! short cut
453#if defined key_lim3
454      INTEGER, DIMENSION(3) ::   zdimsz   ! number of elements in each of the 4 dimensions (i.e. i,j,t,ice-cat) for an array
455      INTEGER               ::   zndims   ! number of dimensions in an array (i.e. 3 = wo ice cat; 4 = w ice cat)
456      INTEGER               ::   inum,id1 ! local integer
457#endif
[3294]458      TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) ::   blf_i         !  array of namelist information structures
459      TYPE(FLD_N) ::   bn_tem, bn_sal, bn_u3d, bn_v3d   !
460      TYPE(FLD_N) ::   bn_ssh, bn_u2d, bn_v2d           ! informations about the fields to be read
[2528]461#if defined key_lim2
[3294]462      TYPE(FLD_N) ::   bn_frld, bn_hicif, bn_hsnif      !
[4292]463#elif defined key_lim3
464      TYPE(FLD_N) ::   bn_a_i, bn_ht_i, bn_ht_s     
[2528]465#endif
[3294]466      NAMELIST/nambdy_dta/ cn_dir, bn_tem, bn_sal, bn_u3d, bn_v3d, bn_ssh, bn_u2d, bn_v2d 
[2528]467#if defined key_lim2
[3294]468      NAMELIST/nambdy_dta/ bn_frld, bn_hicif, bn_hsnif
[4292]469#elif defined key_lim3
470      NAMELIST/nambdy_dta/ bn_a_i, bn_ht_i, bn_ht_s
[3294]471#endif
472      NAMELIST/nambdy_dta/ ln_full_vel
473      !!---------------------------------------------------------------------------
[911]474
[3294]475      IF( nn_timing == 1 ) CALL timing_start('bdy_dta_init')
[911]476
[3651]477      IF(lwp) WRITE(numout,*)
478      IF(lwp) WRITE(numout,*) 'bdy_dta_ini : initialization of data at the open boundaries'
479      IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
480      IF(lwp) WRITE(numout,*) ''
481
[3294]482      ! Set nn_dta
483      DO ib_bdy = 1, nb_bdy
484         nn_dta(ib_bdy) = MAX(  nn_dyn2d_dta(ib_bdy)       &
485                               ,nn_dyn3d_dta(ib_bdy)       &
486                               ,nn_tra_dta(ib_bdy)         &
[4292]487#if ( defined key_lim2 || defined key_lim3 )
488                              ,nn_ice_lim_dta(ib_bdy)    &
[2528]489#endif
[3294]490                              )
491         IF(nn_dta(ib_bdy) .gt. 1) nn_dta(ib_bdy) = 1
492      END DO
[911]493
[3294]494      ! Work out upper bound of how many fields there are to read in and allocate arrays
495      ! ---------------------------------------------------------------------------
496      ALLOCATE( nb_bdy_fld(nb_bdy) )
497      nb_bdy_fld(:) = 0
498      DO ib_bdy = 1, nb_bdy         
[4292]499         IF( cn_dyn2d(ib_bdy) /= 'none' .and. ( nn_dyn2d_dta(ib_bdy) .eq. 1 .or. nn_dyn2d_dta(ib_bdy) .eq. 3 ) ) THEN
[3294]500            nb_bdy_fld(ib_bdy) = nb_bdy_fld(ib_bdy) + 3
501         ENDIF
[4292]502         IF( cn_dyn3d(ib_bdy) /= 'none' .and. nn_dyn3d_dta(ib_bdy) .eq. 1 ) THEN
[3294]503            nb_bdy_fld(ib_bdy) = nb_bdy_fld(ib_bdy) + 2
504         ENDIF
[4292]505         IF( cn_tra(ib_bdy) /= 'none' .and. nn_tra_dta(ib_bdy) .eq. 1  ) THEN
[3294]506            nb_bdy_fld(ib_bdy) = nb_bdy_fld(ib_bdy) + 2
507         ENDIF
[4292]508#if ( defined key_lim2 || defined key_lim3 )
509         IF( cn_ice_lim(ib_bdy) /= 'none' .and. nn_ice_lim_dta(ib_bdy) .eq. 1  ) THEN
[3294]510            nb_bdy_fld(ib_bdy) = nb_bdy_fld(ib_bdy) + 3
511         ENDIF
512#endif               
[3651]513         IF(lwp) WRITE(numout,*) 'Maximum number of files to open =',nb_bdy_fld(ib_bdy)
[3294]514      ENDDO           
[2528]515
[3294]516      nb_bdy_fld_sum = SUM( nb_bdy_fld )
[911]517
[3294]518      ALLOCATE( bf(nb_bdy_fld_sum), STAT=ierror )
519      IF( ierror > 0 ) THEN   
520         CALL ctl_stop( 'bdy_dta: unable to allocate bf structure' )   ;   RETURN 
[911]521      ENDIF
[3294]522      ALLOCATE( blf_i(nb_bdy_fld_sum), STAT=ierror )
523      IF( ierror > 0 ) THEN   
524         CALL ctl_stop( 'bdy_dta: unable to allocate blf_i structure' )   ;   RETURN 
[1125]525      ENDIF
[3294]526      ALLOCATE( nbmap_ptr(nb_bdy_fld_sum), STAT=ierror )
527      IF( ierror > 0 ) THEN   
528         CALL ctl_stop( 'bdy_dta: unable to allocate nbmap_ptr structure' )   ;   RETURN 
529      ENDIF
530      ALLOCATE( ilen1(nb_bdy_fld_sum), ilen3(nb_bdy_fld_sum) ) 
531      ALLOCATE( ibdy(nb_bdy_fld_sum) ) 
532      ALLOCATE( igrid(nb_bdy_fld_sum) ) 
[911]533
[3294]534      ! Read namelists
535      ! --------------
[4147]536      REWIND(numnam_ref)
537      REWIND(numnam_cfg)
[3294]538      jfld = 0 
539      DO ib_bdy = 1, nb_bdy         
540         IF( nn_dta(ib_bdy) .eq. 1 ) THEN
[4147]541            READ  ( numnam_ref, nambdy_dta, IOSTAT = ios, ERR = 901)
542901         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_dta in reference namelist', lwp )
[911]543
[4147]544            READ  ( numnam_cfg, nambdy_dta, IOSTAT = ios, ERR = 902 )
545902         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_dta in configuration namelist', lwp )
546            WRITE ( numond, nambdy_dta )
547
[3294]548            cn_dir_array(ib_bdy) = cn_dir
549            ln_full_vel_array(ib_bdy) = ln_full_vel
[911]550
[3294]551            nblen => idx_bdy(ib_bdy)%nblen
552            nblenrim => idx_bdy(ib_bdy)%nblenrim
[4292]553            dta => dta_bdy(ib_bdy)
554            dta%nread(2) = 0
[911]555
[3294]556            ! Only read in necessary fields for this set.
557            ! Important that barotropic variables come first.
[4292]558            IF( nn_dyn2d_dta(ib_bdy) .eq. 1 .or. nn_dyn2d_dta(ib_bdy) .eq. 3 ) THEN
[911]559
[4292]560               IF( dta%ll_ssh ) THEN
561                  if(lwp) write(numout,*) '++++++ reading in ssh field'
[3294]562                  jfld = jfld + 1
563                  blf_i(jfld) = bn_ssh
564                  ibdy(jfld) = ib_bdy
565                  igrid(jfld) = 1
[3651]566                  ilen1(jfld) = nblen(igrid(jfld))
[3294]567                  ilen3(jfld) = 1
[4292]568                  dta%nread(2) = dta%nread(2) + 1
[3294]569               ENDIF
[911]570
[4292]571               IF( dta%ll_u2d .and. .not. ln_full_vel_array(ib_bdy) ) THEN
572                  if(lwp) write(numout,*) '++++++ reading in u2d field'
[3294]573                  jfld = jfld + 1
574                  blf_i(jfld) = bn_u2d
575                  ibdy(jfld) = ib_bdy
576                  igrid(jfld) = 2
[3651]577                  ilen1(jfld) = nblen(igrid(jfld))
[3294]578                  ilen3(jfld) = 1
[4292]579                  dta%nread(2) = dta%nread(2) + 1
580               ENDIF
[911]581
[4292]582               IF( dta%ll_v2d .and. .not. ln_full_vel_array(ib_bdy) ) THEN
583                  if(lwp) write(numout,*) '++++++ reading in v2d field'
[3294]584                  jfld = jfld + 1
585                  blf_i(jfld) = bn_v2d
586                  ibdy(jfld) = ib_bdy
587                  igrid(jfld) = 3
[3651]588                  ilen1(jfld) = nblen(igrid(jfld))
[3294]589                  ilen3(jfld) = 1
[4292]590                  dta%nread(2) = dta%nread(2) + 1
[3294]591               ENDIF
[911]592
[3294]593            ENDIF
[1125]594
[4292]595            ! read 3D velocities if baroclinic velocities require OR if
596            ! barotropic velocities required and ln_full_vel set to .true.
597            IF( nn_dyn3d_dta(ib_bdy) .eq. 1 .or. &
598           &  ( ln_full_vel_array(ib_bdy) .and. ( nn_dyn2d_dta(ib_bdy) .eq. 1 .or. nn_dyn2d_dta(ib_bdy) .eq. 3 ) ) ) THEN
[911]599
[4292]600               IF( dta%ll_u3d .or. ( ln_full_vel_array(ib_bdy) .and. dta%ll_u2d ) ) THEN
601                  if(lwp) write(numout,*) '++++++ reading in u3d field'
602                  jfld = jfld + 1
603                  blf_i(jfld) = bn_u3d
604                  ibdy(jfld) = ib_bdy
605                  igrid(jfld) = 2
606                  ilen1(jfld) = nblen(igrid(jfld))
607                  ilen3(jfld) = jpk
608                  IF( ln_full_vel_array(ib_bdy) .and. dta%ll_u2d ) dta%nread(2) = dta%nread(2) + 1
609               ENDIF
[911]610
[4292]611               IF( dta%ll_v3d .or. ( ln_full_vel_array(ib_bdy) .and. dta%ll_v2d ) ) THEN
612                  if(lwp) write(numout,*) '++++++ reading in v3d field'
613                  jfld = jfld + 1
614                  blf_i(jfld) = bn_v3d
615                  ibdy(jfld) = ib_bdy
616                  igrid(jfld) = 3
617                  ilen1(jfld) = nblen(igrid(jfld))
618                  ilen3(jfld) = jpk
619                  IF( ln_full_vel_array(ib_bdy) .and. dta%ll_v2d ) dta%nread(2) = dta%nread(2) + 1
620               ENDIF
[911]621
[3294]622            ENDIF
[911]623
[3294]624            ! temperature and salinity
[4292]625            IF( nn_tra_dta(ib_bdy) .eq. 1 ) THEN
[911]626
[4292]627               IF( dta%ll_tem ) THEN
628                  if(lwp) write(numout,*) '++++++ reading in tem field'
629                  jfld = jfld + 1
630                  blf_i(jfld) = bn_tem
631                  ibdy(jfld) = ib_bdy
632                  igrid(jfld) = 1
633                  ilen1(jfld) = nblen(igrid(jfld))
634                  ilen3(jfld) = jpk
635               ENDIF
[911]636
[4292]637               IF( dta%ll_sal ) THEN
638                  if(lwp) write(numout,*) '++++++ reading in sal field'
639                  jfld = jfld + 1
640                  blf_i(jfld) = bn_sal
641                  ibdy(jfld) = ib_bdy
642                  igrid(jfld) = 1
643                  ilen1(jfld) = nblen(igrid(jfld))
644                  ilen3(jfld) = jpk
645               ENDIF
[911]646
647            ENDIF
648
[3294]649#if defined key_lim2
650            ! sea ice
[4292]651            IF( nn_ice_lim2_dta(ib_bdy) .eq. 1 ) THEN
[911]652
[4292]653               IF( dta%ll_frld ) THEN
654                  jfld = jfld + 1
655                  blf_i(jfld) = bn_frld
656                  ibdy(jfld) = ib_bdy
657                  igrid(jfld) = 1
658                  ilen1(jfld) = nblen(igrid(jfld))
659                  ilen3(jfld) = 1
660               ENDIF
[911]661
[4292]662               IF( dta%ll_hicif ) THEN
663                  jfld = jfld + 1
664                  blf_i(jfld) = bn_hicif
665                  ibdy(jfld) = ib_bdy
666                  igrid(jfld) = 1
667                  ilen1(jfld) = nblen(igrid(jfld))
668                  ilen3(jfld) = 1
669               ENDIF
[911]670
[4292]671               IF( dta%ll_hsnif ) THEN
672                  jfld = jfld + 1
673                  blf_i(jfld) = bn_hsnif
674                  ibdy(jfld) = ib_bdy
675                  igrid(jfld) = 1
676                  ilen1(jfld) = nblen(igrid(jfld))
677                  ilen3(jfld) = 1
678               ENDIF
[911]679
[3294]680            ENDIF
[4292]681#elif defined key_lim3
682            ! sea ice
683            IF( nn_ice_lim_dta(ib_bdy) .eq. 1 ) THEN
684
685               ! Test for types of ice input (lim2 or lim3)
686               CALL iom_open ( bn_a_i%clname, inum )
687               id1 = iom_varid ( inum, bn_a_i%clvar, kdimsz=zdimsz, kndims=zndims, ldstop = .FALSE. )
688               CALL iom_close ( inum )
689               !CALL fld_clopn ( bn_a_i, nyear, nmonth, nday, ldstop=.TRUE. )
690               !CALL iom_open ( bn_a_i %clname, inum )
691               !id1 = iom_varid ( bn_a_i%num, bn_a_i%clvar, kdimsz=zdimsz, kndims=zndims, ldstop = .FALSE. )
692                IF ( zndims == 4 ) THEN
693                 ll_bdylim3 = .TRUE.   ! lim3 input
694               ELSE
695                 ll_bdylim3 = .FALSE.  ! lim2 input     
696               ENDIF
697               ! End test
698
699               IF( dta%ll_a_i ) THEN
700                  jfld = jfld + 1
701                  blf_i(jfld) = bn_a_i
702                  ibdy(jfld) = ib_bdy
703                  igrid(jfld) = 1
704                  ilen1(jfld) = nblen(igrid(jfld))
705                  IF ( ll_bdylim3 ) THEN ; ilen3(jfld)=jpl ; ELSE ; ilen3(jfld)=1 ; ENDIF
706               ENDIF
707
708               IF( dta%ll_ht_i ) THEN
709                  jfld = jfld + 1
710                  blf_i(jfld) = bn_ht_i
711                  ibdy(jfld) = ib_bdy
712                  igrid(jfld) = 1
713                  ilen1(jfld) = nblen(igrid(jfld))
714                  IF ( ll_bdylim3 ) THEN ; ilen3(jfld)=jpl ; ELSE ; ilen3(jfld)=1 ; ENDIF
715               ENDIF
716
717               IF( dta%ll_ht_s ) THEN
718                  jfld = jfld + 1
719                   blf_i(jfld) = bn_ht_s
720                  ibdy(jfld) = ib_bdy
721                  igrid(jfld) = 1
722                  ilen1(jfld) = nblen(igrid(jfld))
723                  IF ( ll_bdylim3 ) THEN ; ilen3(jfld)=jpl ; ELSE ; ilen3(jfld)=1 ; ENDIF
724               ENDIF
725
[3294]726#endif
727            ! Recalculate field counts
728            !-------------------------
729            IF( ib_bdy .eq. 1 ) THEN
[4148]730               nb_bdy_fld_sum = 0
[3294]731               nb_bdy_fld(ib_bdy) = jfld
732               nb_bdy_fld_sum     = jfld             
733            ELSE
734               nb_bdy_fld(ib_bdy) = jfld - nb_bdy_fld_sum
735               nb_bdy_fld_sum = nb_bdy_fld_sum + nb_bdy_fld(ib_bdy)
736            ENDIF
[911]737
[4292]738            dta%nread(1) = nb_bdy_fld(ib_bdy)
739
[3294]740         ENDIF ! nn_dta .eq. 1
741      ENDDO ! ib_bdy
[911]742
[3294]743      DO jfld = 1, nb_bdy_fld_sum
744         ALLOCATE( bf(jfld)%fnow(ilen1(jfld),1,ilen3(jfld)) )
745         IF( blf_i(jfld)%ln_tint ) ALLOCATE( bf(jfld)%fdta(ilen1(jfld),1,ilen3(jfld),2) )
746         nbmap_ptr(jfld)%ptr => idx_bdy(ibdy(jfld))%nbmap(:,igrid(jfld))
747      ENDDO
[911]748
[3294]749      ! fill bf with blf_i and control print
750      !-------------------------------------
751      jstart = 1
752      DO ib_bdy = 1, nb_bdy
[3651]753         jend = nb_bdy_fld(ib_bdy) 
[3294]754         CALL fld_fill( bf(jstart:jend), blf_i(jstart:jend), cn_dir_array(ib_bdy), 'bdy_dta',   &
755         &              'open boundary conditions', 'nambdy_dta' )
756         jstart = jend + 1
757      ENDDO
[911]758
[3294]759      ! Initialise local boundary data arrays
760      ! nn_xxx_dta=0 : allocate space - will be filled from initial conditions later
761      ! nn_xxx_dta=1 : point to "fnow" arrays
762      !-------------------------------------
[911]763
[3294]764      jfld = 0
765      DO ib_bdy=1, nb_bdy
[911]766
[3294]767         nblen => idx_bdy(ib_bdy)%nblen
[4292]768         dta => dta_bdy(ib_bdy)
[911]769
[4292]770         if(lwp) then
771            write(numout,*) '++++++ dta%ll_ssh = ',dta%ll_ssh
772            write(numout,*) '++++++ dta%ll_u2d = ',dta%ll_u2d
773            write(numout,*) '++++++ dta%ll_v2d = ',dta%ll_v2d
774            write(numout,*) '++++++ dta%ll_u3d = ',dta%ll_u3d
775            write(numout,*) '++++++ dta%ll_v3d = ',dta%ll_v3d
776            write(numout,*) '++++++ dta%ll_tem = ',dta%ll_tem
777            write(numout,*) '++++++ dta%ll_sal = ',dta%ll_sal
778         endif
779
780         IF ( nn_dyn2d_dta(ib_bdy) .eq. 0 .or. nn_dyn2d_dta(ib_bdy) .eq. 2 ) THEN
781            if(lwp) write(numout,*) '++++++ dta%ssh/u2d/u3d allocated space'
782            IF( dta%ll_ssh ) ALLOCATE( dta%ssh(nblen(1)) )
783            IF( dta%ll_u2d ) ALLOCATE( dta%u2d(nblen(2)) )
784            IF( dta%ll_v2d ) ALLOCATE( dta%v2d(nblen(3)) )
785         ENDIF
786         IF ( nn_dyn2d_dta(ib_bdy) .eq. 1 .or. nn_dyn2d_dta(ib_bdy) .eq. 3 ) THEN
787            IF( dta%ll_ssh ) THEN
788               if(lwp) write(numout,*) '++++++ dta%ssh pointing to fnow'
789               jfld = jfld + 1
790               dta%ssh => bf(jfld)%fnow(:,1,1)
791            ENDIF
792            IF ( dta%ll_u2d ) THEN
793               IF ( ln_full_vel_array(ib_bdy) ) THEN
794                  if(lwp) write(numout,*) '++++++ dta%u2d allocated space'
795                  ALLOCATE( dta%u2d(nblen(2)) )
796               ELSE
797                  if(lwp) write(numout,*) '++++++ dta%u2d pointing to fnow'
[3651]798                  jfld = jfld + 1
[4292]799                  dta%u2d => bf(jfld)%fnow(:,1,1)
800               ENDIF
801            ENDIF
802            IF ( dta%ll_v2d ) THEN
803               IF ( ln_full_vel_array(ib_bdy) ) THEN
804                  if(lwp) write(numout,*) '++++++ dta%v2d allocated space'
805                  ALLOCATE( dta%v2d(nblen(3)) )
[3294]806               ELSE
[4292]807                  if(lwp) write(numout,*) '++++++ dta%v2d pointing to fnow'
[3294]808                  jfld = jfld + 1
[4292]809                  dta%v2d => bf(jfld)%fnow(:,1,1)
[3294]810               ENDIF
811            ENDIF
812         ENDIF
[911]813
[4292]814         IF ( nn_dyn3d_dta(ib_bdy) .eq. 0 ) THEN
815            if(lwp) write(numout,*) '++++++ dta%u3d/v3d allocated space'
816            IF( dta%ll_u3d ) ALLOCATE( dta_bdy(ib_bdy)%u3d(nblen(2),jpk) )
817            IF( dta%ll_v3d ) ALLOCATE( dta_bdy(ib_bdy)%v3d(nblen(3),jpk) )
[3294]818         ENDIF
[4292]819         IF ( nn_dyn3d_dta(ib_bdy) .eq. 1 .or. &
820           &  ( ln_full_vel_array(ib_bdy) .and. ( nn_dyn2d_dta(ib_bdy) .eq. 1 .or. nn_dyn2d_dta(ib_bdy) .eq. 3 ) ) ) THEN
821            IF ( dta%ll_u3d .or. ( ln_full_vel_array(ib_bdy) .and. dta%ll_u2d ) ) THEN
822               if(lwp) write(numout,*) '++++++ dta%u3d pointing to fnow'
823               jfld = jfld + 1
824               dta_bdy(ib_bdy)%u3d => bf(jfld)%fnow(:,1,:)
825            ENDIF
826            IF ( dta%ll_v3d .or. ( ln_full_vel_array(ib_bdy) .and. dta%ll_v2d ) ) THEN
827               if(lwp) write(numout,*) '++++++ dta%v3d pointing to fnow'
828               jfld = jfld + 1
829               dta_bdy(ib_bdy)%v3d => bf(jfld)%fnow(:,1,:)
830            ENDIF
[3294]831         ENDIF
[911]832
[4292]833         IF( nn_tra_dta(ib_bdy) .eq. 0 ) THEN
834            if(lwp) write(numout,*) '++++++ dta%tem/sal allocated space'
835            IF( dta%ll_tem ) ALLOCATE( dta_bdy(ib_bdy)%tem(nblen(1),jpk) )
836            IF( dta%ll_sal ) ALLOCATE( dta_bdy(ib_bdy)%sal(nblen(1),jpk) )
837         ELSE
838            IF( dta%ll_tem ) THEN
839               if(lwp) write(numout,*) '++++++ dta%tem pointing to fnow'
[3294]840               jfld = jfld + 1
841               dta_bdy(ib_bdy)%tem => bf(jfld)%fnow(:,1,:)
[4292]842            ENDIF
843            IF( dta%ll_sal ) THEN
844               if(lwp) write(numout,*) '++++++ dta%sal pointing to fnow'
[3294]845               jfld = jfld + 1
846               dta_bdy(ib_bdy)%sal => bf(jfld)%fnow(:,1,:)
847            ENDIF
848         ENDIF
[911]849
[3294]850#if defined key_lim2
[4292]851         IF (nn_ice_lim(ib_bdy) .gt. 0) THEN
[3294]852            IF( nn_ice_lim2_dta(ib_bdy) .eq. 0 ) THEN
[4292]853               ALLOCATE( dta_bdy(ib_bdy)%frld(nblen(1)) )
854               ALLOCATE( dta_bdy(ib_bdy)%hicif(nblen(1)) )
855               ALLOCATE( dta_bdy(ib_bdy)%hsnif(nblen(1)) )
[3294]856            ELSE
857               jfld = jfld + 1
858               dta_bdy(ib_bdy)%frld  => bf(jfld)%fnow(:,1,1)
859               jfld = jfld + 1
860               dta_bdy(ib_bdy)%hicif => bf(jfld)%fnow(:,1,1)
861               jfld = jfld + 1
862               dta_bdy(ib_bdy)%hsnif => bf(jfld)%fnow(:,1,1)
863            ENDIF
[911]864         ENDIF
[4292]865#elif defined key_lim3
866         IF (nn_ice_lim(ib_bdy) .gt. 0) THEN
867            IF( nn_ice_lim_dta(ib_bdy) .eq. 0 ) THEN
868               ALLOCATE( dta_bdy(ib_bdy)%a_i (nblen(1),jpl) )
869               ALLOCATE( dta_bdy(ib_bdy)%ht_i(nblen(1),jpl) )
870               ALLOCATE( dta_bdy(ib_bdy)%ht_s(nblen(1),jpl) )
871            ELSE
872               IF ( ll_bdylim3 ) THEN ! case input is lim3 type
873                  jfld = jfld + 1
874                  dta_bdy(ib_bdy)%a_i  => bf(jfld)%fnow(:,1,:)
875                  jfld = jfld + 1
876                  dta_bdy(ib_bdy)%ht_i => bf(jfld)%fnow(:,1,:)
877                  jfld = jfld + 1
878                  dta_bdy(ib_bdy)%ht_s => bf(jfld)%fnow(:,1,:)
879               ELSE ! case input is lim2 type
880                  jfld_ai  = jfld + 1
881                  jfld_hti = jfld + 2
882                  jfld_hts = jfld + 3
883                  jfld     = jfld + 3
884                  ALLOCATE( dta_bdy(ib_bdy)%a_i (nblen(1),jpl) )
885                  ALLOCATE( dta_bdy(ib_bdy)%ht_i(nblen(1),jpl) )
886                  ALLOCATE( dta_bdy(ib_bdy)%ht_s(nblen(1),jpl) )
887                  dta_bdy(ib_bdy)%a_i (:,:) = 0.0
888                  dta_bdy(ib_bdy)%ht_i(:,:) = 0.0
889                  dta_bdy(ib_bdy)%ht_s(:,:) = 0.0
890               ENDIF
891
892            ENDIF
893         ENDIF
[3294]894#endif
[911]895
[3294]896      ENDDO ! ib_bdy
[911]897
[3294]898      IF( nn_timing == 1 ) CALL timing_stop('bdy_dta_init')
[911]899
[3294]900      END SUBROUTINE bdy_dta_init
[911]901
902#else
[1125]903   !!----------------------------------------------------------------------
[3294]904   !!   Dummy module                   NO Open Boundary Conditions
[1125]905   !!----------------------------------------------------------------------
[911]906CONTAINS
[3294]907   SUBROUTINE bdy_dta( kt, jit, time_offset ) ! Empty routine
908      INTEGER, INTENT( in )           ::   kt   
909      INTEGER, INTENT( in ), OPTIONAL ::   jit   
910      INTEGER, INTENT( in ), OPTIONAL ::   time_offset
911      WRITE(*,*) 'bdy_dta: You should not have seen this print! error?', kt
912   END SUBROUTINE bdy_dta
913   SUBROUTINE bdy_dta_init()                  ! Empty routine
914      WRITE(*,*) 'bdy_dta_init: You should not have seen this print! error?'
915   END SUBROUTINE bdy_dta_init
[911]916#endif
917
918   !!==============================================================================
919END MODULE bdydta
Note: See TracBrowser for help on using the repository browser.