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.
bdytides.F90 in branches/UKMO/dev_5518_tide_analysis_restart/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/UKMO/dev_5518_tide_analysis_restart/NEMOGCM/NEMO/OPA_SRC/BDY/bdytides.F90 @ 6686

Last change on this file since 6686 was 6686, checked in by deazer, 8 years ago

Added extra option to compute final harmonic analysis or not

File size: 30.7 KB
Line 
1MODULE bdytides
2   !!======================================================================
3   !!                       ***  MODULE  bdytides  ***
4   !! Ocean dynamics:   Tidal forcing at open boundaries
5   !!======================================================================
6   !! History :  2.0  !  2007-01  (D.Storkey)  Original code
7   !!            2.3  !  2008-01  (J.Holt)  Add date correction. Origins POLCOMS v6.3 2007
8   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
9   !!            3.3  !  2010-09  (D.Storkey and E.O'Dea)  bug fixes
10   !!            3.4  !  2012-09  (G. Reffray and J. Chanut) New inputs + mods
11   !!            3.5  !  2013-07  (J. Chanut) Compliant with time splitting changes
12   !!----------------------------------------------------------------------
13#if defined key_bdy
14   !!----------------------------------------------------------------------
15   !!   'key_bdy'     Open Boundary Condition
16   !!----------------------------------------------------------------------
17   !!   PUBLIC
18   !!      bdytide_init     : read of namelist and initialisation of tidal harmonics data
19   !!      tide_update   : calculation of tidal forcing at each timestep
20   !!----------------------------------------------------------------------
21   USE timing          ! Timing
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE iom
25   USE in_out_manager  ! I/O units
26   USE phycst          ! physical constants
27   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
28   USE bdy_par         ! Unstructured boundary parameters
29   USE bdy_oce         ! ocean open boundary conditions
30   USE daymod          ! calendar
31   USE wrk_nemo        ! Memory allocation
32   USE tideini
33!   USE tide_mod       ! Useless ??
34   USE fldread
35   USE dynspg_oce, ONLY: lk_dynspg_ts
36
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   bdytide_init     ! routine called in bdy_init
41   PUBLIC   bdytide_update   ! routine called in bdy_dta
42   PUBLIC   bdy_dta_tides    ! routine called in dyn_spg_ts
43
44   TYPE, PUBLIC ::   TIDES_DATA     !: Storage for external tidal harmonics data
45      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   ssh0       !: Tidal constituents : SSH0 (read in file)
46      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   u0         !: Tidal constituents : U0   (read in file)
47      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   v0         !: Tidal constituents : V0   (read in file)
48      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   ssh        !: Tidal constituents : SSH  (after nodal cor.)
49      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   u          !: Tidal constituents : U    (after nodal cor.)
50      REAL(wp), POINTER, DIMENSION(:,:,:)    ::   v          !: Tidal constituents : V    (after nodal cor.)
51   END TYPE TIDES_DATA
52   INTEGER, PUBLIC, PARAMETER                  ::   jptides_max = 15      !: Max number of tidal contituents
53      LOGICAL, PUBLIC                           ::   ln_harm_ana_store    !: =T Stores data for  harmonic Analysis
54      LOGICAL, PUBLIC                           ::   ln_harm_ana_compute     !: =T  Compute harmonic Analysis
55      LOGICAL, PUBLIC                           ::   ln_harmana_read         !: =T  Decide to do the analysis
56                                                                             !from scratch or continue previous run
57
58!$AGRIF_DO_NOT_TREAT
59   TYPE(TIDES_DATA), PUBLIC, DIMENSION(jp_bdy), TARGET :: tides  !: External tidal harmonics data
60!$AGRIF_END_DO_NOT_TREAT
61   TYPE(OBC_DATA)  , PRIVATE, DIMENSION(jp_bdy) :: dta_bdy_s  !: bdy external data (slow component)
62
63   !!----------------------------------------------------------------------
64   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
65   !! $Id$
66   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
67   !!----------------------------------------------------------------------
68CONTAINS
69
70   SUBROUTINE bdytide_init
71      !!----------------------------------------------------------------------
72      !!                    ***  SUBROUTINE bdytide_init  ***
73      !!                     
74      !! ** Purpose : - Read in namelist for tides and initialise external
75      !!                tidal harmonics data
76      !!
77      !!----------------------------------------------------------------------
78      !! namelist variables
79      !!-------------------
80      CHARACTER(len=80)                         ::   filtide             !: Filename root for tidal input files
81      LOGICAL                                   ::   ln_bdytide_2ddta    !: If true, read 2d harmonic data
82      LOGICAL                                   ::   ln_bdytide_conj     !: If true, assume complex conjugate tidal data
83      !!
84      INTEGER                                   ::   ib_bdy, itide, ib   !: dummy loop indices
85      INTEGER                                   ::   ii, ij              !: dummy loop indices
86      INTEGER                                   ::   inum, igrd
87      INTEGER, DIMENSION(3)                     ::   ilen0       !: length of boundary data (from OBC arrays)
88      INTEGER, POINTER, DIMENSION(:)            ::   nblen, nblenrim     ! short cuts
89      INTEGER                                   ::   ios                 ! Local integer output status for namelist read
90      CHARACTER(len=80)                         ::   clfile              !: full file name for tidal input file
91      REAL(wp),ALLOCATABLE, DIMENSION(:,:,:)    ::   dta_read            !: work space to read in tidal harmonics data
92      REAL(wp), POINTER, DIMENSION(:,:)         ::   ztr, zti            !:  "     "    "   "   "   "        "      "
93      !!
94      TYPE(TIDES_DATA),  POINTER                ::   td                  !: local short cut   
95      TYPE(MAP_POINTER), DIMENSION(jpbgrd)      ::   ibmap_ptr           !: array of pointers to nbmap
96      !!
97      NAMELIST/nambdy_tide/filtide, ln_bdytide_2ddta, ln_bdytide_conj, ln_harm_ana_store, ln_harm_ana_compute, ln_harmana_read
98      !!----------------------------------------------------------------------
99
100      IF( nn_timing == 1 ) CALL timing_start('bdytide_init')
101
102      IF (nb_bdy>0) THEN
103         IF(lwp) WRITE(numout,*)
104         IF(lwp) WRITE(numout,*) 'bdytide_init : initialization of tidal harmonic forcing at open boundaries'
105         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
106      ENDIF
107
108      REWIND(numnam_cfg)
109
110      DO ib_bdy = 1, nb_bdy
111         IF( nn_dyn2d_dta(ib_bdy) .ge. 2 ) THEN
112
113            td => tides(ib_bdy)
114            nblen => idx_bdy(ib_bdy)%nblen
115            nblenrim => idx_bdy(ib_bdy)%nblenrim
116
117            ! Namelist nambdy_tide : tidal harmonic forcing at open boundaries
118            filtide(:) = ''
119
120            ! Don't REWIND here - may need to read more than one of these namelists.
121            READ  ( numnam_ref, nambdy_tide, IOSTAT = ios, ERR = 901)
122901         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_tide in reference namelist', lwp )
123            READ  ( numnam_cfg, nambdy_tide, IOSTAT = ios, ERR = 902 )
124902         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_tide in configuration namelist', lwp )
125            IF(lwm) WRITE ( numond, nambdy_tide )
126            !                                               ! Parameter control and print
127            IF(lwp) WRITE(numout,*) '  '
128            IF(lwp) WRITE(numout,*) '          Namelist nambdy_tide : tidal harmonic forcing at open boundaries'
129            IF(lwp) WRITE(numout,*) '             read tidal data in 2d files: ', ln_bdytide_2ddta
130            IF(lwp) WRITE(numout,*) '             assume complex conjugate   : ', ln_bdytide_conj
131            IF(lwp) WRITE(numout,*) '             Number of tidal components to read: ', nb_harmo
132            IF(lwp) WRITE(numout,*) '             Use PCOMS harmonic ananalysis or not: ', ln_harm_ana_store
133            IF(lwp) WRITE(numout,*) '             Compute Final  harmonic ananalysis or not: ', ln_harm_ana_compute
134            IF(lwp) WRITE(numout,*) '             Read in previous days harmonic data or start afresh: ', ln_harmana_read
135            IF(lwp) THEN
136                    WRITE(numout,*) '             Tidal components: ' 
137               DO itide = 1, nb_harmo
138                  WRITE(numout,*)  '                 ', Wave(ntide(itide))%cname_tide 
139               END DO
140            ENDIF
141            IF(lwp) WRITE(numout,*) ' '
142
143            ! Allocate space for tidal harmonics data - get size from OBC data arrays
144            ! -----------------------------------------------------------------------
145
146            ! JC: If FRS scheme is used, we assume that tidal is needed over the whole
147            ! relaxation area     
148            IF( cn_dyn2d(ib_bdy) == 'frs' ) THEN
149               ilen0(:)=nblen(:)
150            ELSE
151               ilen0(:)=nblenrim(:)
152            ENDIF
153
154            ALLOCATE( td%ssh0( ilen0(1), nb_harmo, 2 ) )
155            ALLOCATE( td%ssh ( ilen0(1), nb_harmo, 2 ) )
156
157            ALLOCATE( td%u0( ilen0(2), nb_harmo, 2 ) )
158            ALLOCATE( td%u ( ilen0(2), nb_harmo, 2 ) )
159
160            ALLOCATE( td%v0( ilen0(3), nb_harmo, 2 ) )
161            ALLOCATE( td%v ( ilen0(3), nb_harmo, 2 ) )
162
163            td%ssh0(:,:,:) = 0._wp
164            td%ssh (:,:,:) = 0._wp
165            td%u0  (:,:,:) = 0._wp
166            td%u   (:,:,:) = 0._wp
167            td%v0  (:,:,:) = 0._wp
168            td%v   (:,:,:) = 0._wp
169
170            IF (ln_bdytide_2ddta) THEN
171               ! It is assumed that each data file contains all complex harmonic amplitudes
172               ! given on the data domain (ie global, jpidta x jpjdta)
173               !
174               CALL wrk_alloc( jpi, jpj, zti, ztr )
175               !
176               ! SSH fields
177               clfile = TRIM(filtide)//'_grid_T.nc'
178               CALL iom_open (clfile , inum ) 
179               igrd = 1                       ! Everything is at T-points here
180               DO itide = 1, nb_harmo
181                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_z1', ztr(:,:) )
182                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_z2', zti(:,:) ) 
183                  DO ib = 1, ilen0(igrd)
184                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
185                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
186                     td%ssh0(ib,itide,1) = ztr(ii,ij)
187                     td%ssh0(ib,itide,2) = zti(ii,ij)
188                  END DO
189               END DO
190               CALL iom_close( inum )
191               !
192               ! U fields
193               clfile = TRIM(filtide)//'_grid_U.nc'
194               CALL iom_open (clfile , inum ) 
195               igrd = 2                       ! Everything is at U-points here
196               DO itide = 1, nb_harmo
197                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_u1', ztr(:,:) )
198                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_u2', zti(:,:) )
199                  DO ib = 1, ilen0(igrd)
200                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
201                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
202                     td%u0(ib,itide,1) = ztr(ii,ij)
203                     td%u0(ib,itide,2) = zti(ii,ij)
204                  END DO
205               END DO
206               CALL iom_close( inum )
207               !
208               ! V fields
209               clfile = TRIM(filtide)//'_grid_V.nc'
210               CALL iom_open (clfile , inum ) 
211               igrd = 3                       ! Everything is at V-points here
212               DO itide = 1, nb_harmo
213                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_v1', ztr(:,:) )
214                  CALL iom_get  ( inum, jpdom_data, TRIM(Wave(ntide(itide))%cname_tide)//'_v2', zti(:,:) )
215                  DO ib = 1, ilen0(igrd)
216                     ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
217                     ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
218                     td%v0(ib,itide,1) = ztr(ii,ij)
219                     td%v0(ib,itide,2) = zti(ii,ij)
220                  END DO
221               END DO 
222               CALL iom_close( inum )
223               !
224               CALL wrk_dealloc( jpi, jpj, ztr, zti ) 
225               !
226            ELSE           
227               !
228               ! Read tidal data only on bdy segments
229               !
230               ALLOCATE( dta_read( MAXVAL(ilen0(1:3)), 1, 1 ) )
231               !
232               ! Set map structure
233               ibmap_ptr(1)%ptr => idx_bdy(ib_bdy)%nbmap(:,1)
234               ibmap_ptr(1)%ll_unstruc = ln_coords_file(ib_bdy)
235               ibmap_ptr(2)%ptr => idx_bdy(ib_bdy)%nbmap(:,2)
236               ibmap_ptr(2)%ll_unstruc = ln_coords_file(ib_bdy)
237               ibmap_ptr(3)%ptr => idx_bdy(ib_bdy)%nbmap(:,3)
238               ibmap_ptr(3)%ll_unstruc = ln_coords_file(ib_bdy)
239
240               ! Open files and read in tidal forcing data
241               ! -----------------------------------------
242
243               DO itide = 1, nb_harmo
244                  !                                                              ! SSH fields
245                  clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_T.nc'
246                  CALL iom_open( clfile, inum )
247                  CALL fld_map( inum, 'z1' , dta_read(1:ilen0(1),1:1,1:1) , 1,  ibmap_ptr(1) )
248                  td%ssh0(:,itide,1) = dta_read(1:ilen0(1),1,1)
249                  CALL fld_map( inum, 'z2' , dta_read(1:ilen0(1),1:1,1:1) , 1,  ibmap_ptr(1) )
250                  td%ssh0(:,itide,2) = dta_read(1:ilen0(1),1,1)
251                  CALL iom_close( inum )
252                  !                                                              ! U fields
253                  clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_U.nc'
254                  CALL iom_open( clfile, inum )
255                  CALL fld_map( inum, 'u1' , dta_read(1:ilen0(2),1:1,1:1) , 1, ibmap_ptr(2) )
256                  td%u0(:,itide,1) = dta_read(1:ilen0(2),1,1)
257                  CALL fld_map( inum, 'u2' , dta_read(1:ilen0(2),1:1,1:1) , 1, ibmap_ptr(2) )
258                  td%u0(:,itide,2) = dta_read(1:ilen0(2),1,1)
259                  CALL iom_close( inum )
260                  !                                                              ! V fields
261                  clfile = TRIM(filtide)//TRIM(Wave(ntide(itide))%cname_tide)//'_grid_V.nc'
262                  CALL iom_open( clfile, inum )
263                  CALL fld_map( inum, 'v1' , dta_read(1:ilen0(3),1:1,1:1) , 1, ibmap_ptr(3) )
264                  td%v0(:,itide,1) = dta_read(1:ilen0(3),1,1)
265                  CALL fld_map( inum, 'v2' , dta_read(1:ilen0(3),1:1,1:1) , 1, ibmap_ptr(3) )
266                  td%v0(:,itide,2) = dta_read(1:ilen0(3),1,1)
267                  CALL iom_close( inum )
268                  !
269               END DO ! end loop on tidal components
270               !
271               DEALLOCATE( dta_read )
272            ENDIF ! ln_bdytide_2ddta=.true.
273            !
274            IF ( ln_bdytide_conj ) THEN ! assume complex conjugate in data files
275               td%ssh0(:,:,2) = - td%ssh0(:,:,2)
276               td%u0  (:,:,2) = - td%u0  (:,:,2)
277               td%v0  (:,:,2) = - td%v0  (:,:,2)
278            ENDIF
279            !
280            IF ( lk_dynspg_ts ) THEN ! Allocate arrays to save slowly varying boundary data during
281                                     ! time splitting integration
282               ALLOCATE( dta_bdy_s(ib_bdy)%ssh ( ilen0(1) ) )
283               ALLOCATE( dta_bdy_s(ib_bdy)%u2d ( ilen0(2) ) )
284               ALLOCATE( dta_bdy_s(ib_bdy)%v2d ( ilen0(3) ) )
285               dta_bdy_s(ib_bdy)%ssh(:) = 0.e0
286               dta_bdy_s(ib_bdy)%u2d(:) = 0.e0
287               dta_bdy_s(ib_bdy)%v2d(:) = 0.e0
288            ENDIF
289            !
290         ENDIF ! nn_dyn2d_dta(ib_bdy) .ge. 2
291         !
292      END DO ! loop on ib_bdy
293
294      IF( nn_timing == 1 ) CALL timing_stop('bdytide_init')
295
296   END SUBROUTINE bdytide_init
297
298   SUBROUTINE bdytide_update ( kt, idx, dta, td, jit, time_offset )
299      !!----------------------------------------------------------------------
300      !!                 ***  SUBROUTINE bdytide_update  ***
301      !!               
302      !! ** Purpose : - Add tidal forcing to ssh, u2d and v2d OBC data arrays.
303      !!               
304      !!----------------------------------------------------------------------
305      INTEGER, INTENT( in )            ::   kt          ! Main timestep counter
306      TYPE(OBC_INDEX), INTENT( in )    ::   idx         ! OBC indices
307      TYPE(OBC_DATA),  INTENT(inout)   ::   dta         ! OBC external data
308      TYPE(TIDES_DATA),INTENT( inout ) ::   td          ! tidal harmonics data
309      INTEGER,INTENT(in),OPTIONAL      ::   jit         ! Barotropic timestep counter (for timesplitting option)
310      INTEGER,INTENT( in ), OPTIONAL   ::   time_offset ! time offset in units of timesteps. NB. if jit
311                                                        ! is present then units = subcycle timesteps.
312                                                        ! time_offset = 0  => get data at "now"    time level
313                                                        ! time_offset = -1 => get data at "before" time level
314                                                        ! time_offset = +1 => get data at "after"  time level
315                                                        ! etc.
316      !!
317      INTEGER, DIMENSION(3)            ::   ilen0       !: length of boundary data (from OBC arrays)
318      INTEGER                          :: itide, igrd, ib   ! dummy loop indices
319      INTEGER                          :: time_add          ! time offset in units of timesteps
320      REAL(wp)                         :: z_arg, z_sarg, zflag, zramp     
321      REAL(wp), DIMENSION(jpmax_harmo) :: z_sist, z_cost
322      !!----------------------------------------------------------------------
323
324      IF( nn_timing == 1 ) CALL timing_start('bdytide_update')
325
326      ilen0(1) =  SIZE(td%ssh(:,1,1))
327      ilen0(2) =  SIZE(td%u(:,1,1))
328      ilen0(3) =  SIZE(td%v(:,1,1))
329
330      zflag=1
331      IF ( PRESENT(jit) ) THEN
332        IF ( jit /= 1 ) zflag=0
333      ENDIF
334
335      IF ( nsec_day == NINT(0.5_wp * rdttra(1)) .AND. zflag==1 ) THEN
336        !
337        kt_tide = kt
338        !
339        IF(lwp) THEN
340           WRITE(numout,*)
341           WRITE(numout,*) 'bdytide_update : (re)Initialization of the tidal bdy forcing at kt=',kt
342           WRITE(numout,*) '~~~~~~~~~~~~~~ '
343        ENDIF
344        !
345        CALL tide_init_elevation ( idx, td )
346        CALL tide_init_velocities( idx, td )
347        !
348      ENDIF
349
350      time_add = 0
351      IF( PRESENT(time_offset) ) THEN
352         time_add = time_offset
353      ENDIF
354         
355      IF( PRESENT(jit) ) THEN 
356         z_arg = ((kt-kt_tide) * rdt + (jit+0.5_wp*(time_add-1)) * rdt / REAL(nn_baro,wp) )
357      ELSE                             
358         z_arg = ((kt-kt_tide)+time_add) * rdt
359      ENDIF
360
361      ! Linear ramp on tidal component at open boundaries
362      zramp = 1._wp
363      IF (ln_tide_ramp) zramp = MIN(MAX( (z_arg + (kt_tide-nit000)*rdt)/(rdttideramp*rday),0._wp),1._wp)
364
365      DO itide = 1, nb_harmo
366         z_sarg = z_arg * omega_tide(itide)
367         z_cost(itide) = COS( z_sarg )
368         z_sist(itide) = SIN( z_sarg )
369      END DO
370
371      DO itide = 1, nb_harmo
372         igrd=1                              ! SSH on tracer grid
373         DO ib = 1, ilen0(igrd)
374            dta%ssh(ib) = dta%ssh(ib) + zramp*(td%ssh(ib,itide,1)*z_cost(itide) + td%ssh(ib,itide,2)*z_sist(itide))
375         END DO
376         igrd=2                              ! U grid
377         DO ib = 1, ilen0(igrd)
378            dta%u2d(ib) = dta%u2d(ib) + zramp*(td%u  (ib,itide,1)*z_cost(itide) + td%u  (ib,itide,2)*z_sist(itide))
379         END DO
380         igrd=3                              ! V grid
381         DO ib = 1, ilen0(igrd) 
382            dta%v2d(ib) = dta%v2d(ib) + zramp*(td%v  (ib,itide,1)*z_cost(itide) + td%v  (ib,itide,2)*z_sist(itide))
383         END DO
384      END DO
385      !
386      IF( nn_timing == 1 ) CALL timing_stop('bdytide_update')
387      !
388   END SUBROUTINE bdytide_update
389
390   SUBROUTINE bdy_dta_tides( kt, kit, time_offset )
391      !!----------------------------------------------------------------------
392      !!                 ***  SUBROUTINE bdy_dta_tides  ***
393      !!               
394      !! ** Purpose : - Add tidal forcing to ssh, u2d and v2d OBC data arrays.
395      !!               
396      !!----------------------------------------------------------------------
397      INTEGER, INTENT( in )            ::   kt          ! Main timestep counter
398      INTEGER, INTENT( in ),OPTIONAL   ::   kit         ! Barotropic timestep counter (for timesplitting option)
399      INTEGER, INTENT( in ),OPTIONAL   ::   time_offset ! time offset in units of timesteps. NB. if kit
400                                                        ! is present then units = subcycle timesteps.
401                                                        ! time_offset = 0  => get data at "now"    time level
402                                                        ! time_offset = -1 => get data at "before" time level
403                                                        ! time_offset = +1 => get data at "after"  time level
404                                                        ! etc.
405      !!
406      LOGICAL  :: lk_first_btstp  ! =.TRUE. if time splitting and first barotropic step
407      INTEGER,          DIMENSION(jpbgrd) :: ilen0 
408      INTEGER, DIMENSION(1:jpbgrd) :: nblen, nblenrim  ! short cuts
409      INTEGER  :: itide, ib_bdy, ib, igrd                     ! loop indices
410      INTEGER  :: time_add                                    ! time offset in units of timesteps
411      REAL(wp) :: z_arg, z_sarg, zramp, zoff, z_cost, z_sist     
412      !!----------------------------------------------------------------------
413
414      IF( nn_timing == 1 ) CALL timing_start('bdy_dta_tides')
415
416      lk_first_btstp=.TRUE.
417      IF ( PRESENT(kit).AND.( kit /= 1 ) ) THEN ; lk_first_btstp=.FALSE. ; ENDIF
418
419      time_add = 0
420      IF( PRESENT(time_offset) ) THEN
421         time_add = time_offset
422      ENDIF
423     
424      ! Absolute time from model initialization:   
425      IF( PRESENT(kit) ) THEN 
426         z_arg = ( kt + (kit+0.5_wp*(time_add-1)) / REAL(nn_baro,wp) ) * rdt
427      ELSE                             
428         z_arg = ( kt + time_add ) * rdt
429      ENDIF
430
431      ! Linear ramp on tidal component at open boundaries
432      zramp = 1.
433      IF (ln_tide_ramp) zramp = MIN(MAX( (z_arg - nit000*rdt)/(rdttideramp*rday),0.),1.)
434
435      DO ib_bdy = 1,nb_bdy
436
437         IF ( nn_dyn2d_dta(ib_bdy) .ge. 2 ) THEN
438
439            nblen(1:jpbgrd) = idx_bdy(ib_bdy)%nblen(1:jpbgrd)
440            nblenrim(1:jpbgrd) = idx_bdy(ib_bdy)%nblenrim(1:jpbgrd)
441
442            IF( cn_dyn2d(ib_bdy) == 'frs' ) THEN
443               ilen0(:)=nblen(:)
444            ELSE
445               ilen0(:)=nblenrim(:)
446            ENDIF     
447
448            ! We refresh nodal factors every day below
449            ! This should be done somewhere else
450            IF ( nsec_day == NINT(0.5_wp * rdttra(1)) .AND. lk_first_btstp ) THEN
451               !
452               kt_tide = kt               
453               !
454               IF(lwp) THEN
455               WRITE(numout,*)
456               WRITE(numout,*) 'bdy_tide_dta : Refresh nodal factors for tidal open bdy data at kt=',kt
457               WRITE(numout,*) '~~~~~~~~~~~~~~ '
458               ENDIF
459               !
460               CALL tide_init_elevation ( idx=idx_bdy(ib_bdy), td=tides(ib_bdy) )
461               CALL tide_init_velocities( idx=idx_bdy(ib_bdy), td=tides(ib_bdy) )
462               !
463            ENDIF
464            zoff = -kt_tide * rdt ! time offset relative to nodal factor computation time
465            !
466            ! If time splitting, save data at first barotropic iteration
467            IF ( PRESENT(kit) ) THEN
468               IF ( lk_first_btstp ) THEN ! Save slow varying open boundary data:
469                  IF ( dta_bdy(ib_bdy)%ll_ssh ) dta_bdy_s(ib_bdy)%ssh(1:ilen0(1)) = dta_bdy(ib_bdy)%ssh(1:ilen0(1))
470                  IF ( dta_bdy(ib_bdy)%ll_u2d ) dta_bdy_s(ib_bdy)%u2d(1:ilen0(2)) = dta_bdy(ib_bdy)%u2d(1:ilen0(2))
471                  IF ( dta_bdy(ib_bdy)%ll_v2d ) dta_bdy_s(ib_bdy)%v2d(1:ilen0(3)) = dta_bdy(ib_bdy)%v2d(1:ilen0(3))
472
473               ELSE ! Initialize arrays from slow varying open boundary data:           
474                  IF ( dta_bdy(ib_bdy)%ll_ssh ) dta_bdy(ib_bdy)%ssh(1:ilen0(1)) = dta_bdy_s(ib_bdy)%ssh(1:ilen0(1))
475                  IF ( dta_bdy(ib_bdy)%ll_u2d ) dta_bdy(ib_bdy)%u2d(1:ilen0(2)) = dta_bdy_s(ib_bdy)%u2d(1:ilen0(2))
476                  IF ( dta_bdy(ib_bdy)%ll_v2d ) dta_bdy(ib_bdy)%v2d(1:ilen0(3)) = dta_bdy_s(ib_bdy)%v2d(1:ilen0(3))
477               ENDIF
478            ENDIF
479            !
480            ! Update open boundary data arrays:
481            DO itide = 1, nb_harmo
482               !
483               z_sarg = (z_arg + zoff) * omega_tide(itide)
484               z_cost = zramp * COS( z_sarg )
485               z_sist = zramp * SIN( z_sarg )
486               !
487               IF ( dta_bdy(ib_bdy)%ll_ssh ) THEN
488                  igrd=1                              ! SSH on tracer grid
489                  DO ib = 1, ilen0(igrd)
490                     dta_bdy(ib_bdy)%ssh(ib) = dta_bdy(ib_bdy)%ssh(ib) + &
491                        &                      ( tides(ib_bdy)%ssh(ib,itide,1)*z_cost + &
492                        &                        tides(ib_bdy)%ssh(ib,itide,2)*z_sist )
493                  END DO
494               ENDIF
495               !
496               IF ( dta_bdy(ib_bdy)%ll_u2d ) THEN
497                  igrd=2                              ! U grid
498                  DO ib = 1, ilen0(igrd)
499                     dta_bdy(ib_bdy)%u2d(ib) = dta_bdy(ib_bdy)%u2d(ib) + &
500                        &                      ( tides(ib_bdy)%u(ib,itide,1)*z_cost + &
501                        &                        tides(ib_bdy)%u(ib,itide,2)*z_sist )
502                  END DO
503               ENDIF
504               !
505               IF ( dta_bdy(ib_bdy)%ll_v2d ) THEN
506                  igrd=3                              ! V grid
507                  DO ib = 1, ilen0(igrd) 
508                     dta_bdy(ib_bdy)%v2d(ib) = dta_bdy(ib_bdy)%v2d(ib) + &
509                        &                      ( tides(ib_bdy)%v(ib,itide,1)*z_cost + &
510                        &                        tides(ib_bdy)%v(ib,itide,2)*z_sist )
511                  END DO
512               ENDIF
513            END DO             
514         END IF
515      END DO
516      !
517      IF( nn_timing == 1 ) CALL timing_stop('bdy_dta_tides')
518      !
519   END SUBROUTINE bdy_dta_tides
520
521   SUBROUTINE tide_init_elevation( idx, td )
522      !!----------------------------------------------------------------------
523      !!                 ***  ROUTINE tide_init_elevation  ***
524      !!----------------------------------------------------------------------
525      TYPE(OBC_INDEX), INTENT( in )      ::   idx     ! OBC indices
526      TYPE(TIDES_DATA),INTENT( inout )   ::   td      ! tidal harmonics data
527      !! * Local declarations
528      INTEGER, DIMENSION(1)            ::   ilen0       !: length of boundary data (from OBC arrays)
529      REAL(wp),ALLOCATABLE, DIMENSION(:) ::   mod_tide, phi_tide
530      INTEGER                            ::   itide, igrd, ib      ! dummy loop indices
531
532      igrd=1   
533                              ! SSH on tracer grid.
534   
535      ilen0(1) =  SIZE(td%ssh0(:,1,1))
536
537      ALLOCATE(mod_tide(ilen0(igrd)),phi_tide(ilen0(igrd)))
538
539      DO itide = 1, nb_harmo
540         DO ib = 1, ilen0(igrd)
541            mod_tide(ib)=SQRT(td%ssh0(ib,itide,1)**2.+td%ssh0(ib,itide,2)**2.)
542            phi_tide(ib)=ATAN2(-td%ssh0(ib,itide,2),td%ssh0(ib,itide,1))
543         END DO
544         DO ib = 1 , ilen0(igrd)
545            mod_tide(ib)=mod_tide(ib)*ftide(itide)
546            phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide)
547         ENDDO
548         DO ib = 1 , ilen0(igrd)
549            td%ssh(ib,itide,1)= mod_tide(ib)*COS(phi_tide(ib))
550            td%ssh(ib,itide,2)=-mod_tide(ib)*SIN(phi_tide(ib))
551         ENDDO
552      END DO
553
554      DEALLOCATE(mod_tide,phi_tide)
555
556   END SUBROUTINE tide_init_elevation
557
558   SUBROUTINE tide_init_velocities( idx, td )
559      !!----------------------------------------------------------------------
560      !!                 ***  ROUTINE tide_init_elevation  ***
561      !!----------------------------------------------------------------------
562      TYPE(OBC_INDEX), INTENT( in )      ::   idx     ! OBC indices
563      TYPE(TIDES_DATA),INTENT( inout )      ::   td      ! tidal harmonics data
564      !! * Local declarations
565      INTEGER, DIMENSION(3)            ::   ilen0       !: length of boundary data (from OBC arrays)
566      REAL(wp),ALLOCATABLE, DIMENSION(:) ::   mod_tide, phi_tide
567      INTEGER                            ::   itide, igrd, ib      ! dummy loop indices
568
569      ilen0(2) =  SIZE(td%u0(:,1,1))
570      ilen0(3) =  SIZE(td%v0(:,1,1))
571
572      igrd=2                                 ! U grid.
573
574      ALLOCATE(mod_tide(ilen0(igrd)),phi_tide(ilen0(igrd)))
575
576      DO itide = 1, nb_harmo
577         DO ib = 1, ilen0(igrd)
578            mod_tide(ib)=SQRT(td%u0(ib,itide,1)**2.+td%u0(ib,itide,2)**2.)
579            phi_tide(ib)=ATAN2(-td%u0(ib,itide,2),td%u0(ib,itide,1))
580         END DO
581         DO ib = 1, ilen0(igrd)
582            mod_tide(ib)=mod_tide(ib)*ftide(itide)
583            phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide)
584         ENDDO
585         DO ib = 1, ilen0(igrd)
586            td%u(ib,itide,1)= mod_tide(ib)*COS(phi_tide(ib))
587            td%u(ib,itide,2)=-mod_tide(ib)*SIN(phi_tide(ib))
588         ENDDO
589      END DO
590
591      DEALLOCATE(mod_tide,phi_tide)
592
593      igrd=3                                 ! V grid.
594
595      ALLOCATE(mod_tide(ilen0(igrd)),phi_tide(ilen0(igrd)))
596
597      DO itide = 1, nb_harmo
598         DO ib = 1, ilen0(igrd)
599            mod_tide(ib)=SQRT(td%v0(ib,itide,1)**2.+td%v0(ib,itide,2)**2.)
600            phi_tide(ib)=ATAN2(-td%v0(ib,itide,2),td%v0(ib,itide,1))
601         END DO
602         DO ib = 1, ilen0(igrd)
603            mod_tide(ib)=mod_tide(ib)*ftide(itide)
604            phi_tide(ib)=phi_tide(ib)+v0tide(itide)+utide(itide)
605         ENDDO
606         DO ib = 1, ilen0(igrd)
607            td%v(ib,itide,1)= mod_tide(ib)*COS(phi_tide(ib))
608            td%v(ib,itide,2)=-mod_tide(ib)*SIN(phi_tide(ib))
609         ENDDO
610      END DO
611
612      DEALLOCATE(mod_tide,phi_tide)
613
614  END SUBROUTINE tide_init_velocities
615#else
616   !!----------------------------------------------------------------------
617   !!   Dummy module         NO Unstruct Open Boundary Conditions for tides
618   !!----------------------------------------------------------------------
619CONTAINS
620   SUBROUTINE bdytide_init             ! Empty routine
621      WRITE(*,*) 'bdytide_init: You should not have seen this print! error?'
622   END SUBROUTINE bdytide_init
623   SUBROUTINE bdytide_update( kt, jit )   ! Empty routine
624      WRITE(*,*) 'bdytide_update: You should not have seen this print! error?', kt, jit
625   END SUBROUTINE bdytide_update
626   SUBROUTINE bdy_dta_tides( kt, kit, time_offset )     ! Empty routine
627      INTEGER, INTENT( in )            ::   kt          ! Dummy argument empty routine     
628      INTEGER, INTENT( in ),OPTIONAL   ::   kit         ! Dummy argument empty routine
629      INTEGER, INTENT( in ),OPTIONAL   ::   time_offset ! Dummy argument empty routine
630      WRITE(*,*) 'bdy_dta_tides: You should not have seen this print! error?', kt, jit
631   END SUBROUTINE bdy_dta_tides
632#endif
633
634   !!======================================================================
635END MODULE bdytides
636
Note: See TracBrowser for help on using the repository browser.