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.
tide_mod.F90 in NEMO/branches/2019/dev_r10742_ENHANCE-12_SimonM-Tides/src/OCE/TDE – NEMO

source: NEMO/branches/2019/dev_r10742_ENHANCE-12_SimonM-Tides/src/OCE/TDE/tide_mod.F90 @ 12056

Last change on this file since 12056 was 12056, checked in by smueller, 4 years ago

Further modifications related to coding style and conventions, addition of citations of equations, and further reduction of the number of real-valued literals used in module tide_mod (ticket #2194)
This changeset affects results produced by the AMM12 and SPITZ12 reference model configurations.

  • Property svn:keywords set to Id
File size: 35.9 KB
Line 
1MODULE tide_mod
2   !!======================================================================
3   !!                       ***  MODULE  tide_mod  ***
4   !! Compute nodal modulations corrections and pulsations
5   !!======================================================================
6   !! History :  1.0  !  2007  (O. Le Galloudec)  Original code
7   !!                 !  2019  (S. Mueller)
8   !!----------------------------------------------------------------------
9   !!
10   !! ** Reference :
11   !!       S58) Schureman, P. (1958): Manual of Harmonic Analysis and
12   !!            Prediction of Tides (Revised (1940) Edition (Reprinted 1958
13   !!            with corrections). Reprinted June 2001). U.S. Department of
14   !!            Commerce, Coast and Geodetic Survey Special Publication
15   !!            No. 98. Washington DC, United States Government Printing
16   !!            Office. 317 pp. DOI: 10.25607/OBP-155.
17   !!----------------------------------------------------------------------
18
19   USE oce, ONLY : sshn       ! sea-surface height
20   USE par_oce                ! ocean parameters
21   USE phycst, ONLY : rpi     ! pi
22   USE daymod, ONLY : ndt05   ! half-length of time step
23   USE in_out_manager         ! I/O units
24   USE iom                    ! xIOs server
25
26   IMPLICIT NONE
27   PRIVATE
28
29   PUBLIC   tide_init
30   PUBLIC   tide_update           ! called by stp
31   PUBLIC   tide_init_harmonics   ! called internally and by module diaharm
32   PUBLIC   upd_tide              ! called in dynspg_... modules
33
34   INTEGER, PUBLIC, PARAMETER ::   jpmax_harmo = 64   !: maximum number of harmonic components
35
36   TYPE ::    tide
37      CHARACTER(LEN=4) ::   cname_tide = ''
38      REAL(wp)         ::   equitide
39      INTEGER          ::   nt, ns, nh, np, np1, shift
40      INTEGER          ::   nksi, nnu0, nnu1, nnu2, R
41      INTEGER          ::   nformula
42   END TYPE tide
43
44   TYPE(tide), DIMENSION(:), POINTER ::   tide_components   !: Array of selected tidal component parameters
45
46   TYPE, PUBLIC ::   tide_harmonic       !:   Oscillation parameters of harmonic tidal components
47      CHARACTER(LEN=4) ::   cname_tide   !    Name of component
48      REAL(wp)         ::   equitide     !    Amplitude of equilibrium tide
49      REAL(wp)         ::   f            !    Node factor
50      REAL(wp)         ::   omega        !    Angular velocity
51      REAL(wp)         ::   v0           !    Initial phase at prime meridian
52      REAL(wp)         ::   u            !    Phase correction
53   END type tide_harmonic
54
55   TYPE(tide_harmonic), PUBLIC, DIMENSION(:), POINTER ::   tide_harmonics   !: Oscillation parameters of selected tidal components
56
57   LOGICAL , PUBLIC ::   ln_tide             !:
58   LOGICAL , PUBLIC ::   ln_tide_pot         !:
59   INTEGER          ::   nn_tide_var         !  Variant of tidal parameter set and tide-potential computation
60   LOGICAL          ::   ln_tide_dia         !  Enable tidal diagnostic output
61   LOGICAL          ::   ln_read_load        !:
62   LOGICAL , PUBLIC ::   ln_scal_load        !:
63   LOGICAL , PUBLIC ::   ln_tide_ramp        !:
64   INTEGER , PUBLIC ::   nb_harmo            !: Number of active tidal components
65   REAL(wp), PUBLIC ::   rn_tide_ramp_dt     !:
66   REAL(wp), PUBLIC ::   rn_scal_load        !:
67   CHARACTER(lc), PUBLIC ::   cn_tide_load   !:
68   REAL(wp)         ::   rn_tide_gamma       ! Tidal tilt factor
69
70   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)   ::   pot_astro            !: tidal potential
71   REAL(wp),         ALLOCATABLE, DIMENSION(:,:)   ::   pot_astro_comp       !  tidal-potential component
72   REAL(wp),         ALLOCATABLE, DIMENSION(:,:,:) ::   amp_pot, phi_pot
73   REAL(wp),         ALLOCATABLE, DIMENSION(:,:,:) ::   amp_load, phi_load
74
75   REAL(wp) ::   rn_tide_ramp_t   !   Elapsed time in seconds
76
77   REAL(wp) ::   sh_T, sh_s, sh_h, sh_p, sh_p1             ! astronomic angles
78   REAL(wp) ::   sh_xi, sh_nu, sh_nuprim, sh_nusec, sh_R   !
79   REAL(wp) ::   sh_I, sh_x1ra, sh_N                       !
80
81   ! Longitudes on 1 Jan 1900, 00h and angular velocities (units of deg and
82   ! deg/h, respectively. The values of these module variables have been copied
83   ! from subroutine astronomic_angle of the version of this module used in
84   ! release version 4.0 of NEMO.
85   REAL(wp) ::   rlon00_N  =  259.1560564_wp               ! Longitude of ascending lunar node
86   REAL(wp) ::   romega_N  = -.0022064139_wp
87   REAL(wp) ::   rlon00_T  =  180.0_wp                     ! Mean solar angle (GMT)
88   REAL(wp) ::   romega_T  =  15.0_wp
89   REAL(wp) ::   rlon00_h  =  280.1895014_wp               ! Mean solar Longitude
90   REAL(wp) ::   romega_h  =  .0410686387_wp
91   REAL(wp) ::   rlon00_s  =  277.0256206_wp               ! Mean lunar Longitude
92   REAL(wp) ::   romega_s  =  .549016532_wp
93   REAL(wp) ::   rlon00_p1 =  281.2208569_wp               ! Longitude of solar perigee
94   REAL(wp) ::   romega_p1 =  .000001961_wp
95   REAL(wp) ::   rlon00_p  =  334.3837214_wp               ! Longitude of lunar perigee
96   REAL(wp) ::   romega_p  =  .004641834_wp
97   ! Values of cos(i)*cos(epsilon), rcice, and sin(incl)*sin(epsilon), rsise,
98   ! where i is the inclination of the orbit of the Moon w.r.t. the ecliptic and
99   ! epsilon the obliquity of the ecliptic on 1 January 1900, 00h. The values of
100   ! these module variables have been copied from subroutine astronomic_angle
101   ! (computation of the cosine of inclination of orbit of Moon to the celestial
102   ! equator) of the version of this module used in release version 4.0 of NEMO.
103   REAL(wp) ::   rcice    = 0.913694997_wp
104   REAL(wp) ::   rsise    = 0.035692561_wp
105   ! Coefficients used to compute sh_xi and sh_nu in subroutine astronomic_angle
106   ! according to two equations given in the explanation of Table 6 of S58
107   REAL(wp) ::   rxinu1, rxinu2
108
109   !!----------------------------------------------------------------------
110   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
111   !! $Id$
112   !! Software governed by the CeCILL license (see ./LICENSE)
113   !!----------------------------------------------------------------------
114CONTAINS
115
116   SUBROUTINE tide_init
117      !!----------------------------------------------------------------------
118      !!                 ***  ROUTINE tide_init  ***
119      !!----------------------------------------------------------------------     
120      INTEGER  :: ji, jk
121      CHARACTER(LEN=4), DIMENSION(jpmax_harmo) :: sn_tide_cnames ! Names of selected tidal components
122      INTEGER  ::   ios                 ! Local integer output status for namelist read
123      !
124      NAMELIST/nam_tide/ln_tide, nn_tide_var, ln_tide_dia, ln_tide_pot, rn_tide_gamma, &
125         &              ln_scal_load, ln_read_load, cn_tide_load,         &
126         &              ln_tide_ramp, rn_scal_load, rn_tide_ramp_dt,      &
127         &              sn_tide_cnames
128      !!----------------------------------------------------------------------
129      !
130      ! Initialise all array elements of sn_tide_cnames, as some of them
131      ! typically do not appear in namelist_ref or namelist_cfg
132      sn_tide_cnames(:) = ''
133      ! Read Namelist nam_tide
134      REWIND( numnam_ref )              ! Namelist nam_tide in reference namelist : Tides
135      READ  ( numnam_ref, nam_tide, IOSTAT = ios, ERR = 901)
136901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_tide in reference namelist', lwp )
137      !
138      REWIND( numnam_cfg )              ! Namelist nam_tide in configuration namelist : Tides
139      READ  ( numnam_cfg, nam_tide, IOSTAT = ios, ERR = 902 )
140902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nam_tide in configuration namelist', lwp )
141      IF(lwm) WRITE ( numond, nam_tide )
142      !
143      IF( ln_tide ) THEN
144         IF (lwp) THEN
145            WRITE(numout,*)
146            WRITE(numout,*) 'tide_init : Initialization of the tidal components'
147            WRITE(numout,*) '~~~~~~~~~ '
148            WRITE(numout,*) '   Namelist nam_tide'
149            WRITE(numout,*) '      Use tidal components                       ln_tide         = ', ln_tide
150            WRITE(numout,*) '         Variant (1: default; 0: legacy option)  nn_tide_var     = ', nn_tide_var
151            WRITE(numout,*) '         Tidal diagnostic output                 ln_tide_dia     = ', ln_tide_dia
152            WRITE(numout,*) '         Apply astronomical potential            ln_tide_pot     = ', ln_tide_pot
153            WRITE(numout,*) '         Tidal tilt factor                       rn_tide_gamma   = ', rn_tide_gamma
154            WRITE(numout,*) '         Use scalar approx. for load potential   ln_scal_load    = ', ln_scal_load
155            WRITE(numout,*) '         Read load potential from file           ln_read_load    = ', ln_read_load
156            WRITE(numout,*) '         Apply ramp on tides at startup          ln_tide_ramp    = ', ln_tide_ramp
157            WRITE(numout,*) '         Fraction of SSH used in scal. approx.   rn_scal_load    = ', rn_scal_load
158            WRITE(numout,*) '         Duration (days) of ramp                 rn_tide_ramp_dt = ', rn_tide_ramp_dt
159         ENDIF
160      ELSE
161         rn_scal_load = 0._wp 
162
163         IF(lwp) WRITE(numout,*)
164         IF(lwp) WRITE(numout,*) 'tide_init : tidal components not used (ln_tide = F)'
165         IF(lwp) WRITE(numout,*) '~~~~~~~~~ '
166         RETURN
167      ENDIF
168      !
169      IF( ln_read_load.AND.(.NOT.ln_tide_pot) ) &
170          &   CALL ctl_stop('ln_read_load requires ln_tide_pot')
171      IF( ln_scal_load.AND.(.NOT.ln_tide_pot) ) &
172          &   CALL ctl_stop('ln_scal_load requires ln_tide_pot')
173      IF( ln_scal_load.AND.ln_read_load ) &
174          &   CALL ctl_stop('Choose between ln_scal_load and ln_read_load')
175      IF( ln_tide_ramp.AND.((nitend-nit000+1)*rdt/rday < rn_tide_ramp_dt) )   &
176         &   CALL ctl_stop('rn_tide_ramp_dt must be lower than run duration')
177      IF( ln_tide_ramp.AND.(rn_tide_ramp_dt<0.) ) &
178         &   CALL ctl_stop('rn_tide_ramp_dt must be positive')
179      !
180      ! Compute coefficients which are used in subroutine astronomic_angle to
181      ! compute sh_xi and sh_nu according to two equations given in the
182      ! explanation of Table 6 of S58
183      rxinu1 = COS( 0.5_wp * ( ABS( ACOS( rcice + rsise ) ) ) ) / COS( 0.5_wp * ( ACOS( rcice - rsise ) ) )
184      rxinu2 = SIN( 0.5_wp * ( ABS( ACOS( rcice + rsise ) ) ) ) / SIN( 0.5_wp * ( ACOS( rcice - rsise ) ) )
185      rxinu1 = 1.01883_wp
186      rxinu2 = 0.64412_wp
187      !
188      ! Initialise array used to store tidal oscillation parameters (frequency,
189      ! amplitude, phase); also retrieve and store array of information about
190      ! selected tidal components
191      CALL tide_init_harmonics(sn_tide_cnames, tide_harmonics, tide_components)
192      !
193      ! Number of active tidal components
194      nb_harmo = size(tide_components)
195      !       
196      ! Ensure that tidal components have been set in namelist_cfg
197      IF( nb_harmo == 0 )   CALL ctl_stop( 'tide_init : No tidal components set in nam_tide' )
198      !
199      IF (.NOT.ln_scal_load ) rn_scal_load = 0._wp
200      !
201      ALLOCATE( amp_pot(jpi,jpj,nb_harmo),                      &
202           &      phi_pot(jpi,jpj,nb_harmo), pot_astro(jpi,jpj)   )
203      IF( ln_tide_dia ) ALLOCATE( pot_astro_comp(jpi,jpj) )
204      IF( ln_read_load ) THEN
205         ALLOCATE( amp_load(jpi,jpj,nb_harmo), phi_load(jpi,jpj,nb_harmo) )
206         CALL tide_init_load
207         amp_pot(:,:,:) = amp_load(:,:,:)
208         phi_pot(:,:,:) = phi_load(:,:,:)
209      ELSE     
210         amp_pot(:,:,:) = 0._wp
211         phi_pot(:,:,:) = 0._wp   
212      ENDIF
213      !
214   END SUBROUTINE tide_init
215
216
217   SUBROUTINE tide_init_components(pcnames, ptide_comp)
218      !!----------------------------------------------------------------------
219      !!                 ***  ROUTINE tide_init_components  ***
220      !!
221      !! Returns pointer to array of variables of type 'tide' that contain
222      !! information about the selected tidal components
223      !! ----------------------------------------------------------------------
224      CHARACTER(LEN=4),              DIMENSION(jpmax_harmo), INTENT(in)  ::   pcnames             ! Names of selected components
225      TYPE(tide),       POINTER,     DIMENSION(:),           INTENT(out) ::   ptide_comp          ! Selected components
226      INTEGER,          ALLOCATABLE, DIMENSION(:)                        ::   icomppos            ! Indices of selected components
227      INTEGER                                                            ::   icomp, jk, jj, ji   ! Miscellaneous integers
228      LOGICAL                                                            ::   llmatch             ! Local variables used for
229      INTEGER                                                            ::   ic1, ic2            !    string comparison
230      TYPE(tide),       POINTER,     DIMENSION(:)                        ::   tide_components     ! All available components
231     
232      ! Populate local array with information about all available tidal
233      ! components
234      !
235      ! Note, here 'tide_components' locally overrides the global module
236      ! variable of the same name to enable the use of the global name in the
237      ! include file that contains the initialisation of elements of array
238      ! 'tide_components'
239      ALLOCATE(tide_components(jpmax_harmo), icomppos(jpmax_harmo))
240      ! Initialise array of indices of the selected componenents
241      icomppos(:) = 0
242      ! Include tidal component parameters for all available components
243      IF (nn_tide_var < 1) THEN
244#define TIDE_VAR_0
245#include "tide.h90"
246#undef TIDE_VAR_0
247      ELSE
248#include "tide.h90"
249      END IF
250      ! Identify the selected components that are availble
251      icomp = 0
252      DO jk = 1, jpmax_harmo
253         IF (TRIM(pcnames(jk)) /= '') THEN
254            DO jj = 1, jpmax_harmo
255               ! Find matches between selected and available constituents
256               ! (ignore capitalisation unless legacy variant has been selected)
257               IF (nn_tide_var < 1) THEN
258                  llmatch = (TRIM(pcnames(jk)) == TRIM(tide_components(jj)%cname_tide))
259               ELSE
260                  llmatch = .TRUE.
261                  ji = MAX(LEN_TRIM(pcnames(jk)), LEN_TRIM(tide_components(jj)%cname_tide))
262                  DO WHILE (llmatch.AND.(ji > 0))
263                     ic1 = IACHAR(pcnames(jk)(ji:ji))
264                     IF ((ic1 >= 97).AND.(ic1 <= 122)) ic1 = ic1 - 32
265                     ic2 = IACHAR(tide_components(jj)%cname_tide(ji:ji))
266                     IF ((ic2 >= 97).AND.(ic2 <= 122)) ic2 = ic2 - 32
267                     llmatch = (ic1 == ic2)
268                     ji = ji - 1
269                  END DO
270               END IF
271               IF (llmatch) THEN
272                  ! Count and record the match
273                  icomp = icomp + 1
274                  icomppos(icomp) = jj
275                  ! Set the capitalisation of the tidal constituent identifier
276                  ! as specified in the namelist
277                  tide_components(jj)%cname_tide = pcnames(jk)
278                  IF (lwp) WRITE(numout, '(10X,"Tidal component #",I2.2,36X,"= ",A4)') icomp, tide_components(jj)%cname_tide
279                  EXIT
280               END IF
281            END DO
282            IF ((lwp).AND.(jj > jpmax_harmo)) WRITE(numout, '(10X,"Tidal component ",A4," is not available!")') pcnames(jk)
283         END IF
284      END DO
285     
286      ! Allocate and populate reduced list of components
287      ALLOCATE(ptide_comp(icomp))
288      DO jk = 1, icomp
289         ptide_comp(jk) = tide_components(icomppos(jk))
290      END DO
291     
292      ! Release local array of available components and list of selected
293      ! components
294      DEALLOCATE(tide_components, icomppos)
295     
296   END SUBROUTINE tide_init_components
297
298
299   SUBROUTINE tide_init_harmonics(pcnames, ptide_harmo, ptide_comp)
300      !!----------------------------------------------------------------------
301      !!                 ***  ROUTINE tide_init_harmonics  ***
302      !!
303      !! Returns pointer to array of variables of type 'tide_harmonics' that
304      !! contain oscillation parameters of the selected harmonic tidal
305      !! components
306      !! ----------------------------------------------------------------------
307      CHARACTER(LEN=4),             DIMENSION(jpmax_harmo), INTENT(in) ::   pcnames     ! Names of selected components
308      TYPE(tide_harmonic), POINTER, DIMENSION(:)                       ::   ptide_harmo ! Oscillation parameters of tidal components
309      TYPE(tide),          POINTER, DIMENSION(:), OPTIONAL             ::   ptide_comp  ! Selected components
310      TYPE(tide),          POINTER, DIMENSION(:)                       ::   ztcomp      ! Selected components
311
312      ! Retrieve information about selected tidal components
313      ! If requested, prepare tidal component array for returning
314      IF (PRESENT(ptide_comp)) THEN
315         CALL tide_init_components(pcnames, ptide_comp)
316         ztcomp => ptide_comp
317      ELSE
318         CALL tide_init_components(pcnames, ztcomp)
319      END IF
320
321      ! Allocate and populate array of oscillation parameters
322      ALLOCATE(ptide_harmo(size(ztcomp)))
323      ptide_harmo(:)%cname_tide = ztcomp(:)%cname_tide
324      ptide_harmo(:)%equitide = ztcomp(:)%equitide
325      CALL tide_harmo(ztcomp, ptide_harmo)
326
327   END SUBROUTINE tide_init_harmonics
328
329
330   SUBROUTINE tide_init_potential
331      !!----------------------------------------------------------------------
332      !!                 ***  ROUTINE tide_init_potential  ***
333      !!
334      !! ** Reference :
335      !!      CT71) Cartwright, D. E. and Tayler, R. J. (1971): New computations of
336      !!            the Tide-generating Potential. Geophys. J. R. astr. Soc. 23,
337      !!            pp. 45-74. DOI: 10.1111/j.1365-246X.1971.tb01803.x
338      !!
339      !!----------------------------------------------------------------------
340      INTEGER  ::   ji, jj, jk   ! dummy loop indices
341      REAL(wp) ::   zcons, ztmp1, ztmp2, zlat, zlon, ztmp, zamp, zcs   ! local scalar
342      !!----------------------------------------------------------------------
343
344      IF( ln_read_load ) THEN
345    amp_pot(:,:,:) = amp_load(:,:,:)
346    phi_pot(:,:,:) = phi_load(:,:,:)
347      ELSE     
348         amp_pot(:,:,:) = 0._wp
349         phi_pot(:,:,:) = 0._wp   
350      ENDIF
351      DO jk = 1, nb_harmo
352         zcons = rn_tide_gamma * tide_components(jk)%equitide * tide_harmonics(jk)%f
353         DO ji = 1, jpi
354            DO jj = 1, jpj
355               ztmp1 =  tide_harmonics(jk)%f * amp_pot(ji,jj,jk) * COS( phi_pot(ji,jj,jk) &
356                  &                                                   + tide_harmonics(jk)%v0 + tide_harmonics(jk)%u )
357               ztmp2 = -tide_harmonics(jk)%f * amp_pot(ji,jj,jk) * SIN( phi_pot(ji,jj,jk) &
358                  &                                                   + tide_harmonics(jk)%v0 + tide_harmonics(jk)%u )
359               zlat = gphit(ji,jj)*rad !! latitude en radian
360               zlon = glamt(ji,jj)*rad !! longitude en radian
361               ztmp = tide_harmonics(jk)%v0 + tide_harmonics(jk)%u + tide_components(jk)%nt * zlon
362               ! le potentiel est composé des effets des astres:
363               SELECT CASE( tide_components(jk)%nt )
364               CASE( 0 )                                                  !   long-periodic tidal constituents (included unless
365                  zcs = zcons * ( 0.5_wp - 1.5_wp * SIN( zlat )**2 )      !   compatibility with original formulation is requested)
366                  IF ( nn_tide_var < 1 ) zcs = 0.0_wp
367               CASE( 1 )                                                  !   diurnal tidal constituents
368                  zcs = zcons * SIN( 2.0_wp*zlat )
369               CASE( 2 )                                                  !   semi-diurnal tidal constituents
370                  zcs = zcons * COS( zlat )**2
371               CASE( 3 )                                                  !   Terdiurnal tidal constituents; the colatitude-dependent
372                  zcs = zcons * COS( zlat )**3                            !   factor is sin(theta)^3 (Table 2 of CT71)
373               CASE DEFAULT                                               !   constituents of higher frequency are not included
374                  zcs = 0.0_wp
375               END SELECT
376               ztmp1 = ztmp1 + zcs * COS( ztmp )
377               ztmp2 = ztmp2 - zcs * SIN( ztmp )
378               zamp = SQRT( ztmp1*ztmp1 + ztmp2*ztmp2 )
379               amp_pot(ji,jj,jk) = zamp
380               phi_pot(ji,jj,jk) = ATAN2( -ztmp2 / MAX( 1.e-10_wp , zamp ) ,   &
381                  &                        ztmp1 / MAX( 1.e-10_wp,  zamp )   )
382            END DO
383         END DO
384      END DO
385      !
386   END SUBROUTINE tide_init_potential
387
388
389   SUBROUTINE tide_init_load
390      !!----------------------------------------------------------------------
391      !!                 ***  ROUTINE tide_init_load  ***
392      !!----------------------------------------------------------------------
393      INTEGER :: inum                 ! Logical unit of input file
394      INTEGER :: ji, jj, itide        ! dummy loop indices
395      REAL(wp), DIMENSION(jpi,jpj) ::   ztr, zti   !: workspace to read in tidal harmonics data
396      !!----------------------------------------------------------------------
397      IF(lwp) THEN
398         WRITE(numout,*)
399         WRITE(numout,*) 'tide_init_load : Initialization of load potential from file'
400         WRITE(numout,*) '~~~~~~~~~~~~~~ '
401      ENDIF
402      !
403      CALL iom_open ( cn_tide_load , inum )
404      !
405      DO itide = 1, nb_harmo
406         CALL iom_get  ( inum, jpdom_data,TRIM(tide_components(itide)%cname_tide)//'_z1', ztr(:,:) )
407         CALL iom_get  ( inum, jpdom_data,TRIM(tide_components(itide)%cname_tide)//'_z2', zti(:,:) )
408         !
409         DO ji=1,jpi
410            DO jj=1,jpj
411               amp_load(ji,jj,itide) =  SQRT( ztr(ji,jj)**2. + zti(ji,jj)**2. )
412               phi_load(ji,jj,itide) = ATAN2(-zti(ji,jj), ztr(ji,jj) )
413            END DO
414         END DO
415         !
416      END DO
417      CALL iom_close( inum )
418      !
419   END SUBROUTINE tide_init_load
420
421
422   SUBROUTINE tide_update( kt )
423      !!----------------------------------------------------------------------
424      !!                 ***  ROUTINE tide_update  ***
425      !!----------------------------------------------------------------------     
426      INTEGER, INTENT( in ) ::   kt     ! ocean time-step
427      INTEGER               ::   jk     ! dummy loop index
428      !!----------------------------------------------------------------------
429     
430      IF( nsec_day == NINT(0.5_wp * rdt) .OR. kt == nit000 ) THEN      ! start a new day
431         !
432         CALL tide_harmo(tide_components, tide_harmonics, ndt05) ! Update oscillation parameters of tidal components for start of current day
433         !
434         !
435         IF(lwp) THEN
436            WRITE(numout,*)
437            WRITE(numout,*) 'tide_update : Update of the components and (re)Init. the potential at kt=', kt
438            WRITE(numout,*) '~~~~~~~~~~~ '
439            DO jk = 1, nb_harmo
440               WRITE(numout,*) tide_harmonics(jk)%cname_tide, tide_harmonics(jk)%u, &
441                  &            tide_harmonics(jk)%f,tide_harmonics(jk)%v0, tide_harmonics(jk)%omega
442            END DO
443         ENDIF
444         !
445         IF( ln_tide_pot )   CALL tide_init_potential
446         !
447         rn_tide_ramp_t = (kt - nit000)*rdt !   Elapsed time in seconds
448      ENDIF
449      !
450   END SUBROUTINE tide_update
451
452
453   SUBROUTINE tide_harmo( ptide_comp, ptide_harmo, psec_day )
454      !
455      TYPE(tide),          DIMENSION(:), POINTER ::   ptide_comp   ! Array of selected tidal component parameters
456      TYPE(tide_harmonic), DIMENSION(:), POINTER ::   ptide_harmo  ! Oscillation parameters of selected tidal components
457      INTEGER, OPTIONAL ::   psec_day                              ! Number of seconds since the start of the current day
458      !
459      IF (PRESENT(psec_day)) THEN
460         CALL astronomic_angle(psec_day)
461      ELSE
462         CALL astronomic_angle(nsec_day)
463      END IF
464      CALL tide_pulse( ptide_comp, ptide_harmo )
465      CALL tide_vuf(   ptide_comp, ptide_harmo )
466      !
467   END SUBROUTINE tide_harmo
468
469
470   SUBROUTINE astronomic_angle(psec_day)
471      !!----------------------------------------------------------------------
472      !!                 ***  ROUTINE astronomic_angle  ***
473      !!                     
474      !! ** Purpose : Compute astronomic angles
475      !!----------------------------------------------------------------------
476      INTEGER  ::   psec_day !   Number of seconds from midnight
477      REAL(wp) ::   zp, zq, zt2, zs2, ztgI2, zP1, ztgn2, zat1, zat2
478      REAL(wp) ::   zqy , zsy, zday, zdj, zhfrac, zt
479      !!----------------------------------------------------------------------
480      !
481      ! Computation of the time from 1 Jan 1900, 00h in years
482      zqy = AINT( (nyear - 1901.0_wp) / 4.0_wp )
483      zsy = nyear - 1900.0_wp
484      !
485      zdj  = dayjul( nyear, nmonth, nday )
486      zday = zdj + zqy - 1.0_wp
487      !
488      zhfrac = psec_day / 3600.0_wp
489      !
490      zt = zsy * 365.0_wp * 24.0_wp + zday * 24.0_wp + zhfrac
491      !
492      ! Longitude of ascending lunar node
493      sh_N  = ( rlon00_N  + romega_N * zt     ) * rad
494      sh_N = MOD( sh_N,  2*rpi )
495      ! Mean solar angle (Greenwhich time)
496      sh_T  = ( rlon00_T  + romega_T * zhfrac ) * rad
497      ! Mean solar Longitude
498      sh_h  = ( rlon00_h  + romega_h * zt     ) * rad
499      sh_h = MOD( sh_h,  2*rpi )
500      ! Mean lunar Longitude
501      sh_s  = ( rlon00_s  + romega_s * zt     ) * rad
502      sh_s = MOD( sh_s,  2*rpi )
503      ! Longitude of solar perigee
504      sh_p1 = ( rlon00_p1 + romega_p1 * zt    ) * rad
505      sh_p1= MOD( sh_p1, 2*rpi )
506      ! Longitude of lunar perigee
507      sh_p  = ( rlon00_p  + romega_p  * zt    ) * rad
508      sh_p = MOD( sh_p,  2*rpi )
509      !
510      ! Inclination of the orbit of the moon w.r.t. the celestial equator, see
511      ! explanation of Table 6 of S58
512      sh_I = ACOS( rcice - rsise * COS( sh_N ) )
513      !
514      ! Computation of sh_xi and sh_nu, see explanation of Table 6 of S58
515      ztgn2 = TAN( sh_N / 2.0_wp )
516      zat1 = ATAN( rxinu1 * ztgn2 )
517      zat2 = ATAN( rxinu2 * ztgn2 )
518      sh_xi = sh_N - zat1 - zat2
519      IF( sh_N > rpi ) sh_xi = sh_xi - 2.0_wp * rpi
520      sh_nu = zat1 - zat2
521      !
522      ! Computation of sh_x1ra, sh_R, sh_nuprim, and sh_nusec used for tidal
523      ! constituents L2, K1, and K2
524      !
525      ! Computation of sh_x1ra and sh_R (Equations 204, 213, and 214 of S58)
526      ztgI2 = tan( sh_I / 2.0_wp )
527      zP1   = sh_p - sh_xi
528      zt2   = ztgI2 * ztgI2
529      sh_x1ra = SQRT( 1.0 - 12.0 * zt2 * COS( 2.0_wp * zP1 ) + 36.0_wp * zt2 * zt2 )
530      zp = SIN( 2.0_wp * zP1 )
531      zq = 1.0_wp / ( 6.0_wp * zt2 ) - COS( 2.0_wp * zP1 )
532      sh_R = ATAN( zp / zq )
533      !
534      ! Computation of sh_nuprim (Equation 224 of S58)
535      zp = SIN( 2.0_wp * sh_I ) * SIN( sh_nu )
536      zq = SIN( 2.0_wp * sh_I ) * COS( sh_nu ) + 0.3347_wp
537      sh_nuprim = ATAN( zp / zq )
538      !
539      ! Computation of sh_nusec  (Equation 232 of S58)
540      zs2 = SIN( sh_I ) * SIN( sh_I )
541      zp  = zs2 * SIN( 2.0_wp * sh_nu )
542      zq  = zs2 * COS( 2.0_wp * sh_nu ) + 0.0727_wp
543      sh_nusec = 0.5_wp * ATAN( zp / zq )
544      !
545   END SUBROUTINE astronomic_angle
546
547
548   SUBROUTINE tide_pulse( ptide_comp, ptide_harmo )
549      !!----------------------------------------------------------------------
550      !!                     ***  ROUTINE tide_pulse  ***
551      !!                     
552      !! ** Purpose : Compute tidal frequencies
553      !!----------------------------------------------------------------------
554      TYPE(tide),          DIMENSION(:), POINTER ::   ptide_comp   ! Array of selected tidal component parameters
555      TYPE(tide_harmonic), DIMENSION(:), POINTER ::   ptide_harmo  ! Oscillation parameters of selected tidal components
556      !
557      INTEGER  ::   jh
558      REAL(wp) ::   zscale
559      !!----------------------------------------------------------------------
560      !
561      zscale =  rad / 3600.0_wp
562      !
563      DO jh = 1, size(ptide_harmo)
564         ptide_harmo(jh)%omega = (  romega_T * ptide_comp( jh )%nT   &
565            &                     + romega_s * ptide_comp( jh )%ns   &
566            &                     + romega_h * ptide_comp( jh )%nh   &
567            &                     + romega_p * ptide_comp( jh )%np   &
568            &                     + romega_p1* ptide_comp( jh )%np1  ) * zscale
569      END DO
570      !
571   END SUBROUTINE tide_pulse
572
573
574   SUBROUTINE tide_vuf( ptide_comp, ptide_harmo )
575      !!----------------------------------------------------------------------
576      !!                     ***  ROUTINE tide_vuf  ***
577      !!                     
578      !! ** Purpose : Compute nodal modulation corrections
579      !!
580      !! ** Outputs : vt: Phase of tidal potential relative to Greenwich (radians)
581      !!              ut: Phase correction u due to nodal motion (radians)
582      !!              ft: Nodal correction factor
583      !!----------------------------------------------------------------------
584      TYPE(tide),          DIMENSION(:), POINTER ::   ptide_comp   ! Array of selected tidal component parameters
585      TYPE(tide_harmonic), DIMENSION(:), POINTER ::   ptide_harmo  ! Oscillation parameters of selected tidal components
586      !
587      INTEGER ::   jh   ! dummy loop index
588      !!----------------------------------------------------------------------
589      !
590      DO jh = 1, size(ptide_harmo)
591         !  Phase of the tidal potential relative to the Greenwhich
592         !  meridian (e.g. the position of the fictuous celestial body). Units are radian:
593         ptide_harmo(jh)%v0 = sh_T * ptide_comp( jh )%nT    &
594            &                   + sh_s * ptide_comp( jh )%ns    &
595            &                   + sh_h * ptide_comp( jh )%nh    &
596            &                   + sh_p * ptide_comp( jh )%np    &
597            &                   + sh_p1* ptide_comp( jh )%np1   &
598            &                   +        ptide_comp( jh )%shift * rad
599         !
600         !  Phase correction u due to nodal motion. Units are radian:
601         ptide_harmo(jh)%u = sh_xi     * ptide_comp( jh )%nksi   &
602            &                  + sh_nu     * ptide_comp( jh )%nnu0   &
603            &                  + sh_nuprim * ptide_comp( jh )%nnu1   &
604            &                  + sh_nusec  * ptide_comp( jh )%nnu2   &
605            &                  + sh_R      * ptide_comp( jh )%R
606
607         !  Nodal correction factor:
608         ptide_harmo(jh)%f = nodal_factort( ptide_comp( jh )%nformula )
609      END DO
610      !
611   END SUBROUTINE tide_vuf
612
613
614   RECURSIVE FUNCTION nodal_factort( kformula ) RESULT( zf )
615      !!----------------------------------------------------------------------
616      !!                  ***  FUNCTION nodal_factort  ***
617      !!
618      !! ** Purpose : Compute amplitude correction factors due to nodal motion
619      !!----------------------------------------------------------------------
620      INTEGER, INTENT(in) ::   kformula
621      !
622      REAL(wp)            ::   zf
623      REAL(wp)            ::   zs, zf1, zf2
624      CHARACTER(LEN=3)    ::   clformula
625      !!----------------------------------------------------------------------
626      !
627      SELECT CASE( kformula )
628      !
629      CASE( 0 )                  ! Formula 0, solar waves
630         zf = 1.0
631         !
632      CASE( 1 )                  ! Formula 1, compound waves (78 x 78)
633         zf=nodal_factort( 78 )
634         zf = zf * zf
635         !
636      CASE ( 4 )                 ! Formula 4,  compound waves (78 x 235)
637         zf1 = nodal_factort( 78 )
638         zf  = nodal_factort(235)
639         zf  = zf1 * zf
640         !
641      CASE( 18 )                 ! Formula 18,  compound waves (78 x 78 x 78 )
642         zf1 = nodal_factort( 78 )
643         zf  = zf1 * zf1 * zf1
644         !
645      CASE( 20 )                 ! Formula 20, compound waves ( 78 x 78 x 78 x 78 )
646         zf1 = nodal_factort( 78 )
647         zf  = zf1 * zf1 * zf1 * zf1
648         !
649      CASE( 73 )                 ! Formula 73 of S58
650         zs = SIN( sh_I )
651         zf = ( 2.0_wp / 3.0_wp - zs * zs ) / 0.5021_wp
652         !
653      CASE( 74 )                 ! Formula 74 of S58
654         zs = SIN(sh_I)
655         zf = zs * zs / 0.1578_wp
656         !
657      CASE( 75 )                 ! Formula 75 of S58
658         zs = COS( sh_I / 2.0_wp )
659         zf = SIN( sh_I ) * zs * zs / 0.3800_wp
660         !
661      CASE( 76 )                 ! Formula 76 of S58
662         zf = SIN( 2.0_wp * sh_I ) / 0.7214_wp
663         !
664      CASE( 78 )                 ! Formula 78 of S58
665         zs = COS( sh_I/2 )
666         zf = zs * zs * zs * zs / 0.9154_wp
667         !
668      CASE( 149 )                ! Formula 149 of S58
669         zs = COS( sh_I/2 )
670         zf = zs * zs * zs * zs * zs * zs / 0.8758_wp
671         !
672      CASE( 215 )                ! Formula 215 of S58 with typo correction (0.9154 instead of 0.9145)
673         zs = COS( sh_I/2 )
674         zf = zs * zs * zs * zs / 0.9154_wp * sh_x1ra
675         !
676      CASE( 227 )                ! Formula 227 of S58
677         zs = SIN( 2.0_wp * sh_I )
678         zf = SQRT( 0.8965_wp * zs * zs + 0.6001_wp * zs * COS( sh_nu ) + 0.1006_wp )
679         !
680      CASE ( 235 )               ! Formula 235 of S58
681         zs = SIN( sh_I )
682         zf = SQRT( 19.0444_wp * zs * zs * zs * zs + 2.7702_wp * zs * zs * cos( 2.0_wp * sh_nu ) + 0.0981_wp )
683         !
684      CASE DEFAULT
685         WRITE( clformula, '(I3)' ) kformula
686         CALL ctl_stop('nodal_factort: formula ' // clformula // ' is not available')
687      END SELECT
688      !
689   END FUNCTION nodal_factort
690
691
692   FUNCTION dayjul( kyr, kmonth, kday )
693      !!----------------------------------------------------------------------
694      !!                   ***  FUNCTION dayjul  ***
695      !!
696      !! Purpose :  compute the Julian day
697      !!----------------------------------------------------------------------
698      INTEGER,INTENT(in)    ::   kyr, kmonth, kday
699      !
700      INTEGER,DIMENSION(12) ::   idayt = (/ 0, 31, 59, 90, 120, 151, 181, 212, 243, 273, 304, 334 /)
701      INTEGER,DIMENSION(12) ::   idays
702      INTEGER               ::   inc, ji, zyq
703      REAL(wp)              ::   dayjul
704      !!----------------------------------------------------------------------
705      !
706      idays(1) = 0
707      idays(2) = 31
708      inc = 0.0_wp
709      zyq = MOD( kyr - 1900 , 4 )
710      IF( zyq == 0 ) inc = 1
711      DO ji = 3, 12
712         idays(ji) = idayt(ji) + inc
713      END DO
714      dayjul = REAL( idays(kmonth) + kday, KIND=wp )
715      !
716   END FUNCTION dayjul
717
718
719   SUBROUTINE upd_tide(pdelta)
720      !!----------------------------------------------------------------------
721      !!                 ***  ROUTINE upd_tide  ***
722      !!
723      !! ** Purpose :   provide at each time step the astronomical potential
724      !!
725      !! ** Method  :   computed from pulsation and amplitude of all tide components
726      !!
727      !! ** Action  :   pot_astro   actronomical potential
728      !!----------------------------------------------------------------------     
729      REAL, INTENT(in)              ::   pdelta      ! Temporal offset in seconds
730      INTEGER                       ::   jk          ! Dummy loop index
731      REAL(wp)                      ::   zt, zramp   ! Local scalars
732      REAL(wp), DIMENSION(nb_harmo) ::   zwt         ! Temporary array
733      !!----------------------------------------------------------------------     
734      !
735      zwt(:) = tide_harmonics(:)%omega * pdelta
736      !
737      IF( ln_tide_ramp ) THEN         ! linear increase if asked
738         zt = rn_tide_ramp_t + pdelta
739         zramp = MIN(  MAX( zt / (rn_tide_ramp_dt*rday) , 0._wp ) , 1._wp  )
740      ENDIF
741      !
742      pot_astro(:,:) = 0._wp          ! update tidal potential (sum of all harmonics)
743      DO jk = 1, nb_harmo
744         IF ( .NOT. ln_tide_dia ) THEN
745            pot_astro(:,:) = pot_astro(:,:) + amp_pot(:,:,jk) * COS( zwt(jk) + phi_pot(:,:,jk) )
746         ELSE
747            pot_astro_comp(:,:) = amp_pot(:,:,jk) * COS( zwt(jk) + phi_pot(:,:,jk) )
748            pot_astro(:,:) = pot_astro(:,:) + pot_astro_comp(:,:)
749            IF ( iom_use( "tide_pot_" // TRIM( tide_harmonics(jk)%cname_tide ) ) ) THEN   ! Output tidal potential (incl. load potential)
750               IF ( ln_tide_ramp ) pot_astro_comp(:,:) = zramp * pot_astro_comp(:,:)
751               CALL iom_put( "tide_pot_" // TRIM( tide_harmonics(jk)%cname_tide ), pot_astro_comp(:,:) )
752            END IF
753         END IF
754      END DO
755      !
756      IF ( ln_tide_ramp ) pot_astro(:,:) = zramp * pot_astro(:,:)
757      !
758      IF( ln_tide_dia ) THEN          ! Output total tidal potential (incl. load potential)
759         IF ( iom_use( "tide_pot" ) ) CALL iom_put( "tide_pot", pot_astro(:,:) + rn_scal_load * sshn(:,:) )
760      END IF
761      !
762   END SUBROUTINE upd_tide
763
764   !!======================================================================
765END MODULE tide_mod
Note: See TracBrowser for help on using the repository browser.