tide_mod.F90 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411
  1. MODULE 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 dom_oce ! ocean space and time domain
  9. USE phycst ! physical constant
  10. USE daymod ! calendar
  11. IMPLICIT NONE
  12. PRIVATE
  13. PUBLIC tide_harmo ! called by tideini and diaharm modules
  14. PUBLIC tide_init_Wave ! called by tideini and diaharm modules
  15. INTEGER, PUBLIC, PARAMETER :: jpmax_harmo = 19 !: maximum number of harmonic
  16. TYPE, PUBLIC :: tide
  17. CHARACTER(LEN=4) :: cname_tide
  18. REAL(wp) :: equitide
  19. INTEGER :: nutide
  20. INTEGER :: nt, ns, nh, np, np1, shift
  21. INTEGER :: nksi, nnu0, nnu1, nnu2, R
  22. INTEGER :: nformula
  23. END TYPE tide
  24. TYPE(tide), PUBLIC, DIMENSION(jpmax_harmo) :: Wave !:
  25. REAL(wp) :: sh_T, sh_s, sh_h, sh_p, sh_p1 ! astronomic angles
  26. REAL(wp) :: sh_xi, sh_nu, sh_nuprim, sh_nusec, sh_R !
  27. REAL(wp) :: sh_I, sh_x1ra, sh_N !
  28. !!----------------------------------------------------------------------
  29. !! NEMO/OPA 3.3 , LOCEAN-IPSL (2010)
  30. !! $Id: tide_mod.F90 2355 2015-05-20 07:11:50Z ufla $
  31. !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
  32. !!----------------------------------------------------------------------
  33. CONTAINS
  34. SUBROUTINE tide_init_Wave
  35. # include "tide.h90"
  36. END SUBROUTINE tide_init_Wave
  37. SUBROUTINE tide_harmo( pomega, pvt, put , pcor, ktide ,kc)
  38. !!----------------------------------------------------------------------
  39. !!----------------------------------------------------------------------
  40. INTEGER , DIMENSION(kc), INTENT(in ) :: ktide ! Indice of tidal constituents
  41. INTEGER , INTENT(in ) :: kc ! Total number of tidal constituents
  42. REAL(wp), DIMENSION(kc), INTENT(out) :: pomega ! pulsation in radians/s
  43. REAL(wp), DIMENSION(kc), INTENT(out) :: pvt, put, pcor !
  44. !!----------------------------------------------------------------------
  45. !
  46. CALL astronomic_angle
  47. CALL tide_pulse( pomega, ktide ,kc )
  48. CALL tide_vuf ( pvt, put, pcor, ktide ,kc )
  49. !
  50. END SUBROUTINE tide_harmo
  51. SUBROUTINE astronomic_angle
  52. !!----------------------------------------------------------------------
  53. !! tj is time elapsed since 1st January 1900, 0 hour, counted in julian
  54. !! century (e.g. time in days divide by 36525)
  55. !!----------------------------------------------------------------------
  56. REAL(wp) :: cosI, p, q, t2, t4, sin2I, s2, tgI2, P1, sh_tgn2, at1, at2
  57. REAL(wp) :: zqy , zsy, zday, zdj, zhfrac
  58. !!----------------------------------------------------------------------
  59. !
  60. zqy = AINT( (nyear-1901.)/4. )
  61. zsy = nyear - 1900.
  62. !
  63. zdj = dayjul( nyear, nmonth, nday )
  64. zday = zdj + zqy - 1.
  65. !
  66. zhfrac = nsec_day / 3600.
  67. !
  68. !----------------------------------------------------------------------
  69. ! Sh_n Longitude of ascending lunar node
  70. !----------------------------------------------------------------------
  71. sh_N=(259.1560564-19.328185764*zsy-.0529539336*zday-.0022064139*zhfrac)*rad
  72. !----------------------------------------------------------------------
  73. ! T mean solar angle (Greenwhich time)
  74. !----------------------------------------------------------------------
  75. sh_T=(180.+zhfrac*(360./24.))*rad
  76. !----------------------------------------------------------------------
  77. ! h mean solar Longitude
  78. !----------------------------------------------------------------------
  79. sh_h=(280.1895014-.238724988*zsy+.9856473288*zday+.0410686387*zhfrac)*rad
  80. !----------------------------------------------------------------------
  81. ! s mean lunar Longitude
  82. !----------------------------------------------------------------------
  83. sh_s=(277.0256206+129.38482032*zsy+13.176396768*zday+.549016532*zhfrac)*rad
  84. !----------------------------------------------------------------------
  85. ! p1 Longitude of solar perigee
  86. !----------------------------------------------------------------------
  87. sh_p1=(281.2208569+.01717836*zsy+.000047064*zday+.000001961*zhfrac)*rad
  88. !----------------------------------------------------------------------
  89. ! p Longitude of lunar perigee
  90. !----------------------------------------------------------------------
  91. sh_p=(334.3837214+40.66246584*zsy+.111404016*zday+.004641834*zhfrac)*rad
  92. sh_N = MOD( sh_N ,2*rpi )
  93. sh_s = MOD( sh_s ,2*rpi )
  94. sh_h = MOD( sh_h, 2*rpi )
  95. sh_p = MOD( sh_p, 2*rpi )
  96. sh_p1= MOD( sh_p1,2*rpi )
  97. cosI = 0.913694997 -0.035692561 *cos(sh_N)
  98. sh_I = ACOS( cosI )
  99. sin2I = sin(sh_I)
  100. sh_tgn2 = tan(sh_N/2.0)
  101. at1=atan(1.01883*sh_tgn2)
  102. at2=atan(0.64412*sh_tgn2)
  103. sh_xi=-at1-at2+sh_N
  104. IF( sh_N > rpi ) sh_xi=sh_xi-2.0*rpi
  105. sh_nu = at1 - at2
  106. !----------------------------------------------------------------------
  107. ! For constituents l2 k1 k2
  108. !----------------------------------------------------------------------
  109. tgI2 = tan(sh_I/2.0)
  110. P1 = sh_p-sh_xi
  111. t2 = tgI2*tgI2
  112. t4 = t2*t2
  113. sh_x1ra = sqrt( 1.0-12.0*t2*cos(2.0*P1)+36.0*t4 )
  114. p = sin(2.0*P1)
  115. q = 1.0/(6.0*t2)-cos(2.0*P1)
  116. sh_R = atan(p/q)
  117. p = sin(2.0*sh_I)*sin(sh_nu)
  118. q = sin(2.0*sh_I)*cos(sh_nu)+0.3347
  119. sh_nuprim = atan(p/q)
  120. s2 = sin(sh_I)*sin(sh_I)
  121. p = s2*sin(2.0*sh_nu)
  122. q = s2*cos(2.0*sh_nu)+0.0727
  123. sh_nusec = 0.5*atan(p/q)
  124. !
  125. END SUBROUTINE astronomic_angle
  126. SUBROUTINE tide_pulse( pomega, ktide ,kc )
  127. !!----------------------------------------------------------------------
  128. !! *** ROUTINE tide_pulse ***
  129. !!
  130. !! ** Purpose : Compute tidal frequencies
  131. !!----------------------------------------------------------------------
  132. INTEGER , INTENT(in ) :: kc ! Total number of tidal constituents
  133. INTEGER , DIMENSION(kc), INTENT(in ) :: ktide ! Indice of tidal constituents
  134. REAL(wp), DIMENSION(kc), INTENT(out) :: pomega ! pulsation in radians/s
  135. !
  136. INTEGER :: jh
  137. REAL(wp) :: zscale
  138. REAL(wp) :: zomega_T = 13149000.0_wp
  139. REAL(wp) :: zomega_s = 481267.892_wp
  140. REAL(wp) :: zomega_h = 36000.76892_wp
  141. REAL(wp) :: zomega_p = 4069.0322056_wp
  142. REAL(wp) :: zomega_n = 1934.1423972_wp
  143. REAL(wp) :: zomega_p1= 1.719175_wp
  144. !!----------------------------------------------------------------------
  145. !
  146. zscale = rad / ( 36525._wp * 86400._wp )
  147. !
  148. DO jh = 1, kc
  149. pomega(jh) = ( zomega_T * Wave( ktide(jh) )%nT &
  150. & + zomega_s * Wave( ktide(jh) )%ns &
  151. & + zomega_h * Wave( ktide(jh) )%nh &
  152. & + zomega_p * Wave( ktide(jh) )%np &
  153. & + zomega_p1* Wave( ktide(jh) )%np1 ) * zscale
  154. END DO
  155. !
  156. END SUBROUTINE tide_pulse
  157. SUBROUTINE tide_vuf( pvt, put, pcor, ktide ,kc )
  158. !!----------------------------------------------------------------------
  159. !! *** ROUTINE tide_vuf ***
  160. !!
  161. !! ** Purpose : Compute nodal modulation corrections
  162. !!
  163. !! ** Outputs : vt: Phase of tidal potential relative to Greenwich (radians)
  164. !! ut: Phase correction u due to nodal motion (radians)
  165. !! ft: Nodal correction factor
  166. !!----------------------------------------------------------------------
  167. INTEGER , INTENT(in ) :: kc ! Total number of tidal constituents
  168. INTEGER , DIMENSION(kc), INTENT(in ) :: ktide ! Indice of tidal constituents
  169. REAL(wp), DIMENSION(kc), INTENT(out) :: pvt, put, pcor !
  170. !
  171. INTEGER :: jh ! dummy loop index
  172. !!----------------------------------------------------------------------
  173. !
  174. DO jh = 1, kc
  175. ! Phase of the tidal potential relative to the Greenwhich
  176. ! meridian (e.g. the position of the fictuous celestial body). Units are radian:
  177. pvt(jh) = sh_T * Wave( ktide(jh) )%nT &
  178. & + sh_s * Wave( ktide(jh) )%ns &
  179. & + sh_h * Wave( ktide(jh) )%nh &
  180. & + sh_p * Wave( ktide(jh) )%np &
  181. & + sh_p1* Wave( ktide(jh) )%np1 &
  182. & + Wave( ktide(jh) )%shift * rad
  183. !
  184. ! Phase correction u due to nodal motion. Units are radian:
  185. put(jh) = sh_xi * Wave( ktide(jh) )%nksi &
  186. & + sh_nu * Wave( ktide(jh) )%nnu0 &
  187. & + sh_nuprim * Wave( ktide(jh) )%nnu1 &
  188. & + sh_nusec * Wave( ktide(jh) )%nnu2 &
  189. & + sh_R * Wave( ktide(jh) )%R
  190. ! Nodal correction factor:
  191. pcor(jh) = nodal_factort( Wave( ktide(jh) )%nformula )
  192. END DO
  193. !
  194. END SUBROUTINE tide_vuf
  195. RECURSIVE FUNCTION nodal_factort( kformula ) RESULT( zf )
  196. !!----------------------------------------------------------------------
  197. !!----------------------------------------------------------------------
  198. INTEGER, INTENT(in) :: kformula
  199. !
  200. REAL(wp) :: zf
  201. REAL(wp) :: zs, zf1, zf2
  202. !!----------------------------------------------------------------------
  203. !
  204. SELECT CASE( kformula )
  205. !
  206. CASE( 0 ) !== formule 0, solar waves
  207. zf = 1.0
  208. !
  209. CASE( 1 ) !== formule 1, compound waves (78 x 78)
  210. zf=nodal_factort(78)
  211. zf = zf * zf
  212. !
  213. CASE ( 2 ) !== formule 2, compound waves (78 x 0) === (78)
  214. zf1= nodal_factort(78)
  215. zf = nodal_factort( 0)
  216. zf = zf1 * zf
  217. !
  218. CASE ( 4 ) !== formule 4, compound waves (78 x 235)
  219. zf1 = nodal_factort( 78)
  220. zf = nodal_factort(235)
  221. zf = zf1 * zf
  222. !
  223. CASE ( 5 ) !== formule 5, compound waves (78 *78 x 235)
  224. zf1 = nodal_factort( 78)
  225. zf = nodal_factort(235)
  226. zf = zf * zf1 * zf1
  227. !
  228. CASE ( 6 ) !== formule 6, compound waves (78 *78 x 0)
  229. zf1 = nodal_factort(78)
  230. zf = nodal_factort( 0)
  231. zf = zf * zf1 * zf1
  232. !
  233. CASE( 7 ) !== formule 7, compound waves (75 x 75)
  234. zf = nodal_factort(75)
  235. zf = zf * zf
  236. !
  237. CASE( 8 ) !== formule 8, compound waves (78 x 0 x 235)
  238. zf = nodal_factort( 78)
  239. zf1 = nodal_factort( 0)
  240. zf2 = nodal_factort(235)
  241. zf = zf * zf1 * zf2
  242. !
  243. CASE( 9 ) !== formule 9, compound waves (78 x 0 x 227)
  244. zf = nodal_factort( 78)
  245. zf1 = nodal_factort( 0)
  246. zf2 = nodal_factort(227)
  247. zf = zf * zf1 * zf2
  248. !
  249. CASE( 10 ) !== formule 10, compound waves (78 x 227)
  250. zf = nodal_factort( 78)
  251. zf1 = nodal_factort(227)
  252. zf = zf * zf1
  253. !
  254. CASE( 11 ) !== formule 11, compound waves (75 x 0)
  255. !!gm bug???? zf 2 fois !
  256. zf = nodal_factort(75)
  257. zf = nodal_factort( 0)
  258. zf = zf * zf1
  259. !
  260. CASE( 12 ) !== formule 12, compound waves (78 x 78 x 78 x 0)
  261. zf1 = nodal_factort(78)
  262. zf = nodal_factort( 0)
  263. zf = zf * zf1 * zf1 * zf1
  264. !
  265. CASE( 13 ) !== formule 13, compound waves (78 x 75)
  266. zf1 = nodal_factort(78)
  267. zf = nodal_factort(75)
  268. zf = zf * zf1
  269. !
  270. CASE( 14 ) !== formule 14, compound waves (235 x 0) === (235)
  271. zf = nodal_factort(235)
  272. zf1 = nodal_factort( 0)
  273. zf = zf * zf1
  274. !
  275. CASE( 15 ) !== formule 15, compound waves (235 x 75)
  276. zf = nodal_factort(235)
  277. zf1 = nodal_factort( 75)
  278. zf = zf * zf1
  279. !
  280. CASE( 16 ) !== formule 16, compound waves (78 x 0 x 0) === (78)
  281. zf = nodal_factort(78)
  282. zf1 = nodal_factort( 0)
  283. zf = zf * zf1 * zf1
  284. !
  285. CASE( 17 ) !== formule 17, compound waves (227 x 0)
  286. zf1 = nodal_factort(227)
  287. zf = nodal_factort( 0)
  288. zf = zf * zf1
  289. !
  290. CASE( 18 ) !== formule 18, compound waves (78 x 78 x 78 )
  291. zf1 = nodal_factort(78)
  292. zf = zf1 * zf1 * zf1
  293. !
  294. CASE( 19 ) !== formule 19, compound waves (78 x 0 x 0 x 0) === (78)
  295. !!gm bug2 ==>>> here identical to formule 16, a third multiplication by zf1 is missing
  296. zf = nodal_factort(78)
  297. zf1 = nodal_factort( 0)
  298. zf = zf * zf1 * zf1
  299. !
  300. CASE( 73 ) !== formule 73
  301. zs = sin(sh_I)
  302. zf = (2./3.-zs*zs)/0.5021
  303. !
  304. CASE( 74 ) !== formule 74
  305. zs = sin(sh_I)
  306. zf = zs * zs / 0.1578
  307. !
  308. CASE( 75 ) !== formule 75
  309. zs = cos(sh_I/2)
  310. zf = sin(sh_I) * zs * zs / 0.3800
  311. !
  312. CASE( 76 ) !== formule 76
  313. zf = sin(2*sh_I) / 0.7214
  314. !
  315. CASE( 77 ) !== formule 77
  316. zs = sin(sh_I/2)
  317. zf = sin(sh_I) * zs * zs / 0.0164
  318. !
  319. CASE( 78 ) !== formule 78
  320. zs = cos(sh_I/2)
  321. zf = zs * zs * zs * zs / 0.9154
  322. !
  323. CASE( 79 ) !== formule 79
  324. zs = sin(sh_I)
  325. zf = zs * zs / 0.1565
  326. !
  327. CASE( 144 ) !== formule 144
  328. zs = sin(sh_I/2)
  329. zf = ( 1-10*zs*zs+15*zs*zs*zs*zs ) * cos(sh_I/2) / 0.5873
  330. !
  331. CASE( 149 ) !== formule 149
  332. zs = cos(sh_I/2)
  333. zf = zs*zs*zs*zs*zs*zs / 0.8758
  334. !
  335. CASE( 215 ) !== formule 215
  336. zs = cos(sh_I/2)
  337. zf = zs*zs*zs*zs / 0.9154 * sh_x1ra
  338. !
  339. CASE( 227 ) !== formule 227
  340. zs = sin(2*sh_I)
  341. zf = sqrt( 0.8965*zs*zs+0.6001*zs*cos (sh_nu)+0.1006 )
  342. !
  343. CASE ( 235 ) !== formule 235
  344. zs = sin(sh_I)
  345. zf = sqrt( 19.0444*zs*zs*zs*zs + 2.7702*zs*zs*cos(2*sh_nu) + .0981 )
  346. !
  347. END SELECT
  348. !
  349. END FUNCTION nodal_factort
  350. FUNCTION dayjul( kyr, kmonth, kday )
  351. !!----------------------------------------------------------------------
  352. !! *** THIS ROUTINE COMPUTES THE JULIAN DAY (AS A REAL VARIABLE)
  353. !!----------------------------------------------------------------------
  354. INTEGER,INTENT(in) :: kyr, kmonth, kday
  355. !
  356. INTEGER,DIMENSION(12) :: idayt, idays
  357. INTEGER :: inc, ji
  358. REAL(wp) :: dayjul, zyq
  359. !
  360. DATA idayt/0.,31.,59.,90.,120.,151.,181.,212.,243.,273.,304.,334./
  361. !!----------------------------------------------------------------------
  362. !
  363. idays(1) = 0.
  364. idays(2) = 31.
  365. inc = 0.
  366. zyq = MOD( kyr-1900. , 4. )
  367. IF( zyq == 0.) inc = 1.
  368. DO ji = 3, 12
  369. idays(ji)=idayt(ji)+inc
  370. END DO
  371. dayjul = idays(kmonth) + kday
  372. !
  373. END FUNCTION dayjul
  374. !!======================================================================
  375. END MODULE tide_mod