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 @ 11704

Last change on this file since 11704 was 11704, checked in by smueller, 5 years ago

Expansion of include file tide.h90 to provide a choice of alternative (currently identical) tidal-constituent parameter sets during code preprocessing (ticket #2194)

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