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

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

Simplification of initialisation of tidal-constituent parameters (module tide_mod) (ticket #2194)

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