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

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

Addition of terdiurnal components to the computation of the tidal potential (ticket #2194)

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