1 | MODULE bdytides |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE bdytides *** |
---|
4 | !! Ocean dynamics: Tidal forcing at open boundaries |
---|
5 | !!====================================================================== |
---|
6 | !! History : 2.0 ! 2007-01 (D.Storkey) Original code |
---|
7 | !! 2.3 ! 2008-01 (J.Holt) Add date correction. Origins POLCOMS v6.3 2007 |
---|
8 | !! 3.0 ! 2008-04 (NEMO team) add in the reference version |
---|
9 | !! 3.3 ! 2010-09 (D.Storkey and E.O'Dea) bug fixes |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | #if defined key_bdy |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! 'key_bdy' Open Boundary Condition |
---|
14 | !!---------------------------------------------------------------------- |
---|
15 | !! PUBLIC |
---|
16 | !! tide_init : read of namelist and initialisation of tidal harmonics data |
---|
17 | !! tide_update : calculation of tidal forcing at each timestep |
---|
18 | !! PRIVATE |
---|
19 | !! uvset :\ |
---|
20 | !! vday : | Routines to correct tidal harmonics forcing for |
---|
21 | !! shpen : | start time of integration |
---|
22 | !! ufset : | |
---|
23 | !! vset :/ |
---|
24 | !!---------------------------------------------------------------------- |
---|
25 | USE oce ! ocean dynamics and tracers |
---|
26 | USE dom_oce ! ocean space and time domain |
---|
27 | USE iom |
---|
28 | USE in_out_manager ! I/O units |
---|
29 | USE phycst ! physical constants |
---|
30 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
31 | USE bdy_par ! Unstructured boundary parameters |
---|
32 | USE bdy_oce ! ocean open boundary conditions |
---|
33 | USE daymod ! calendar |
---|
34 | USE fldread, ONLY: fld_map |
---|
35 | |
---|
36 | IMPLICIT NONE |
---|
37 | PRIVATE |
---|
38 | |
---|
39 | PUBLIC tide_init ! routine called in nemo_init |
---|
40 | PUBLIC tide_update ! routine called in bdydyn |
---|
41 | |
---|
42 | TYPE, PUBLIC :: TIDES_DATA !: Storage for external tidal harmonics data |
---|
43 | INTEGER :: ncpt !: Actual number of tidal components |
---|
44 | REAL(wp), POINTER, DIMENSION(:) :: speed !: Phase speed of tidal constituent (deg/hr) |
---|
45 | REAL(wp), POINTER, DIMENSION(:,:,:) :: ssh !: Tidal constituents : SSH |
---|
46 | REAL(wp), POINTER, DIMENSION(:,:,:) :: u !: Tidal constituents : U |
---|
47 | REAL(wp), POINTER, DIMENSION(:,:,:) :: v !: Tidal constituents : V |
---|
48 | END TYPE TIDES_DATA |
---|
49 | |
---|
50 | INTEGER, PUBLIC, PARAMETER :: jptides_max = 15 !: Max number of tidal contituents |
---|
51 | |
---|
52 | TYPE(TIDES_DATA), PUBLIC, DIMENSION(jp_bdy), TARGET :: tides !: External tidal harmonics data |
---|
53 | |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
56 | !! $Id$ |
---|
57 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
58 | !!---------------------------------------------------------------------- |
---|
59 | CONTAINS |
---|
60 | |
---|
61 | SUBROUTINE tide_init |
---|
62 | !!---------------------------------------------------------------------- |
---|
63 | !! *** SUBROUTINE tide_init *** |
---|
64 | !! |
---|
65 | !! ** Purpose : - Read in namelist for tides and initialise external |
---|
66 | !! tidal harmonics data |
---|
67 | !! |
---|
68 | !!---------------------------------------------------------------------- |
---|
69 | !! namelist variables |
---|
70 | !!------------------- |
---|
71 | LOGICAL :: ln_tide_date !: =T correct tide phases and amplitude for model start date |
---|
72 | CHARACTER(len=80) :: filtide !: Filename root for tidal input files |
---|
73 | CHARACTER(len= 4), DIMENSION(jptides_max) :: tide_cpt !: Names of tidal components used. |
---|
74 | REAL(wp), DIMENSION(jptides_max) :: tide_speed !: Phase speed of tidal constituent (deg/hr) |
---|
75 | !! |
---|
76 | INTEGER, DIMENSION(jptides_max) :: nindx !: index of pre-set tidal components |
---|
77 | INTEGER :: ib_bdy, itide, ib !: dummy loop indices |
---|
78 | INTEGER :: inum, igrd |
---|
79 | INTEGER, DIMENSION(3) :: ilen0 !: length of boundary data (from OBC arrays) |
---|
80 | CHARACTER(len=80) :: clfile !: full file name for tidal input file |
---|
81 | REAL(wp) :: z_arg, z_atde, z_btde, z1t, z2t |
---|
82 | REAL(wp),DIMENSION(jptides_max) :: z_vplu, z_ftc |
---|
83 | REAL(wp),ALLOCATABLE, DIMENSION(:,:,:) :: dta_read !: work space to read in tidal harmonics data |
---|
84 | !! |
---|
85 | TYPE(TIDES_DATA), POINTER :: td !: local short cut |
---|
86 | !! |
---|
87 | NAMELIST/nambdy_tide/ln_tide_date, filtide, tide_cpt, tide_speed |
---|
88 | !!---------------------------------------------------------------------- |
---|
89 | |
---|
90 | IF(lwp) WRITE(numout,*) |
---|
91 | IF(lwp) WRITE(numout,*) 'tide_init : initialization of tidal harmonic forcing at open boundaries' |
---|
92 | IF(lwp) WRITE(numout,*) '~~~~~~~~~' |
---|
93 | |
---|
94 | REWIND(numnam) |
---|
95 | DO ib_bdy = 1, nb_bdy |
---|
96 | IF( nn_dyn2d_dta(ib_bdy) .ge. 2 ) THEN |
---|
97 | |
---|
98 | td => tides(ib_bdy) |
---|
99 | |
---|
100 | ! Namelist nambdy_tide : tidal harmonic forcing at open boundaries |
---|
101 | ln_tide_date = .false. |
---|
102 | filtide(:) = '' |
---|
103 | tide_speed(:) = 0.0 |
---|
104 | tide_cpt(:) = '' |
---|
105 | |
---|
106 | ! Don't REWIND here - may need to read more than one of these namelists. |
---|
107 | READ ( numnam, nambdy_tide ) |
---|
108 | ! ! Count number of components specified |
---|
109 | td%ncpt = 0 |
---|
110 | DO itide = 1, jptides_max |
---|
111 | IF( tide_cpt(itide) /= '' ) THEN |
---|
112 | td%ncpt = td%ncpt + 1 |
---|
113 | ENDIF |
---|
114 | END DO |
---|
115 | |
---|
116 | ! Fill in phase speeds from namelist |
---|
117 | ALLOCATE( td%speed(td%ncpt) ) |
---|
118 | td%speed = tide_speed(1:td%ncpt) |
---|
119 | |
---|
120 | ! Find constituents in standard list |
---|
121 | DO itide = 1, td%ncpt |
---|
122 | nindx(itide) = 0 |
---|
123 | IF( TRIM( tide_cpt(itide) ) == 'Q1' ) nindx(itide) = 1 |
---|
124 | IF( TRIM( tide_cpt(itide) ) == 'O1' ) nindx(itide) = 2 |
---|
125 | IF( TRIM( tide_cpt(itide) ) == 'P1' ) nindx(itide) = 3 |
---|
126 | IF( TRIM( tide_cpt(itide) ) == 'S1' ) nindx(itide) = 4 |
---|
127 | IF( TRIM( tide_cpt(itide) ) == 'K1' ) nindx(itide) = 5 |
---|
128 | IF( TRIM( tide_cpt(itide) ) == '2N2' ) nindx(itide) = 6 |
---|
129 | IF( TRIM( tide_cpt(itide) ) == 'MU2' ) nindx(itide) = 7 |
---|
130 | IF( TRIM( tide_cpt(itide) ) == 'N2' ) nindx(itide) = 8 |
---|
131 | IF( TRIM( tide_cpt(itide) ) == 'NU2' ) nindx(itide) = 9 |
---|
132 | IF( TRIM( tide_cpt(itide) ) == 'M2' ) nindx(itide) = 10 |
---|
133 | IF( TRIM( tide_cpt(itide) ) == 'L2' ) nindx(itide) = 11 |
---|
134 | IF( TRIM( tide_cpt(itide) ) == 'T2' ) nindx(itide) = 12 |
---|
135 | IF( TRIM( tide_cpt(itide) ) == 'S2' ) nindx(itide) = 13 |
---|
136 | IF( TRIM( tide_cpt(itide) ) == 'K2' ) nindx(itide) = 14 |
---|
137 | IF( TRIM( tide_cpt(itide) ) == 'M4' ) nindx(itide) = 15 |
---|
138 | IF( nindx(itide) == 0 .AND. lwp ) THEN |
---|
139 | WRITE(ctmp1,*) 'constitunent', itide,':', tide_cpt(itide), 'not in standard list' |
---|
140 | CALL ctl_warn( ctmp1 ) |
---|
141 | ENDIF |
---|
142 | END DO |
---|
143 | ! ! Parameter control and print |
---|
144 | IF( td%ncpt < 1 ) THEN |
---|
145 | CALL ctl_stop( ' Did not find any tidal components in namelist nambdy_tide' ) |
---|
146 | ELSE |
---|
147 | IF(lwp) WRITE(numout,*) ' Namelist nambdy_tide : tidal harmonic forcing at open boundaries' |
---|
148 | IF(lwp) WRITE(numout,*) ' tidal components specified ', td%ncpt |
---|
149 | IF(lwp) WRITE(numout,*) ' ', tide_cpt(1:td%ncpt) |
---|
150 | IF(lwp) WRITE(numout,*) ' associated phase speeds (deg/hr) : ' |
---|
151 | IF(lwp) WRITE(numout,*) ' ', tide_speed(1:td%ncpt) |
---|
152 | ENDIF |
---|
153 | |
---|
154 | |
---|
155 | ! Allocate space for tidal harmonics data - |
---|
156 | ! get size from OBC data arrays |
---|
157 | ! --------------------------------------- |
---|
158 | |
---|
159 | ilen0(1) = SIZE( dta_bdy(ib_bdy)%ssh ) |
---|
160 | ALLOCATE( td%ssh( ilen0(1), td%ncpt, 2 ) ) |
---|
161 | |
---|
162 | ilen0(2) = SIZE( dta_bdy(ib_bdy)%u2d ) |
---|
163 | ALLOCATE( td%u( ilen0(2), td%ncpt, 2 ) ) |
---|
164 | |
---|
165 | ilen0(3) = SIZE( dta_bdy(ib_bdy)%v2d ) |
---|
166 | ALLOCATE( td%v( ilen0(3), td%ncpt, 2 ) ) |
---|
167 | |
---|
168 | ALLOCATE( dta_read( MAXVAL(ilen0), 1, 1 ) ) |
---|
169 | |
---|
170 | |
---|
171 | ! Open files and read in tidal forcing data |
---|
172 | ! ----------------------------------------- |
---|
173 | |
---|
174 | DO itide = 1, td%ncpt |
---|
175 | ! ! SSH fields |
---|
176 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_T.nc' |
---|
177 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
178 | CALL iom_open( clfile, inum ) |
---|
179 | CALL fld_map( inum, 'z1' , dta_read(1:ilen0(1),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,1) ) |
---|
180 | td%ssh(:,itide,1) = dta_read(1:ilen0(1),1,1) |
---|
181 | CALL fld_map( inum, 'z2' , dta_read(1:ilen0(1),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,1) ) |
---|
182 | td%ssh(:,itide,2) = dta_read(1:ilen0(1),1,1) |
---|
183 | CALL iom_close( inum ) |
---|
184 | ! ! U fields |
---|
185 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_U.nc' |
---|
186 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
187 | CALL iom_open( clfile, inum ) |
---|
188 | CALL fld_map( inum, 'u1' , dta_read(1:ilen0(2),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,2) ) |
---|
189 | td%u(:,itide,1) = dta_read(1:ilen0(2),1,1) |
---|
190 | CALL fld_map( inum, 'u2' , dta_read(1:ilen0(2),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,2) ) |
---|
191 | td%u(:,itide,2) = dta_read(1:ilen0(2),1,1) |
---|
192 | CALL iom_close( inum ) |
---|
193 | ! ! V fields |
---|
194 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_V.nc' |
---|
195 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
196 | CALL iom_open( clfile, inum ) |
---|
197 | CALL fld_map( inum, 'v1' , dta_read(1:ilen0(3),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,3) ) |
---|
198 | td%v(:,itide,1) = dta_read(1:ilen0(3),1,1) |
---|
199 | CALL fld_map( inum, 'v2' , dta_read(1:ilen0(3),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,3) ) |
---|
200 | td%v(:,itide,2) = dta_read(1:ilen0(3),1,1) |
---|
201 | CALL iom_close( inum ) |
---|
202 | ! |
---|
203 | END DO ! end loop on tidal components |
---|
204 | |
---|
205 | IF( ln_tide_date ) THEN ! correct for date factors |
---|
206 | |
---|
207 | !! used nmonth, nyear and nday from daymod.... |
---|
208 | ! Calculate date corrects for 15 standard consituents |
---|
209 | ! This is the initialisation step, so nday, nmonth, nyear are the |
---|
210 | ! initial date/time of the integration. |
---|
211 | print *, nday,nmonth,nyear |
---|
212 | nyear = int(ndate0 / 10000 ) ! initial year |
---|
213 | nmonth = int((ndate0 - nyear * 10000 ) / 100 ) ! initial month |
---|
214 | nday = int(ndate0 - nyear * 10000 - nmonth * 100) |
---|
215 | |
---|
216 | CALL uvset( 0, nday, nmonth, nyear, z_ftc, z_vplu ) |
---|
217 | |
---|
218 | IF(lwp) WRITE(numout,*) 'Correcting tide for date:', nday, nmonth, nyear |
---|
219 | |
---|
220 | DO itide = 1, td%ncpt ! loop on tidal components |
---|
221 | ! |
---|
222 | IF( nindx(itide) /= 0 ) THEN |
---|
223 | !!gm use rpi and rad global variable |
---|
224 | z_arg = 3.14159265d0 * z_vplu(nindx(itide)) / 180.0d0 |
---|
225 | z_atde=z_ftc(nindx(itide))*cos(z_arg) |
---|
226 | z_btde=z_ftc(nindx(itide))*sin(z_arg) |
---|
227 | IF(lwp) WRITE(numout,'(2i5,8f10.6)') itide, nindx(itide), td%speed(itide), & |
---|
228 | & z_ftc(nindx(itide)), z_vplu(nindx(itide)) |
---|
229 | ELSE |
---|
230 | z_atde = 1.0_wp |
---|
231 | z_btde = 0.0_wp |
---|
232 | ENDIF |
---|
233 | ! ! elevation |
---|
234 | igrd = 1 |
---|
235 | DO ib = 1, ilen0(igrd) |
---|
236 | z1t = z_atde * td%ssh(ib,itide,1) + z_btde * td%ssh(ib,itide,2) |
---|
237 | z2t = z_atde * td%ssh(ib,itide,2) - z_btde * td%ssh(ib,itide,1) |
---|
238 | td%ssh(ib,itide,1) = z1t |
---|
239 | td%ssh(ib,itide,2) = z2t |
---|
240 | END DO |
---|
241 | ! ! u |
---|
242 | igrd = 2 |
---|
243 | DO ib = 1, ilen0(igrd) |
---|
244 | z1t = z_atde * td%u(ib,itide,1) + z_btde * td%u(ib,itide,2) |
---|
245 | z2t = z_atde * td%u(ib,itide,2) - z_btde * td%u(ib,itide,1) |
---|
246 | td%u(ib,itide,1) = z1t |
---|
247 | td%u(ib,itide,2) = z2t |
---|
248 | END DO |
---|
249 | ! ! v |
---|
250 | igrd = 3 |
---|
251 | DO ib = 1, ilen0(igrd) |
---|
252 | z1t = z_atde * td%v(ib,itide,1) + z_btde * td%v(ib,itide,2) |
---|
253 | z2t = z_atde * td%v(ib,itide,2) - z_btde * td%v(ib,itide,1) |
---|
254 | td%v(ib,itide,1) = z1t |
---|
255 | td%v(ib,itide,2) = z2t |
---|
256 | END DO |
---|
257 | ! |
---|
258 | END DO ! end loop on tidal components |
---|
259 | ! |
---|
260 | ENDIF ! date correction |
---|
261 | ! |
---|
262 | ENDIF ! nn_dyn2d_dta(ib_bdy) .ge. 2 |
---|
263 | ! |
---|
264 | END DO ! loop on ib_bdy |
---|
265 | |
---|
266 | END SUBROUTINE tide_init |
---|
267 | |
---|
268 | |
---|
269 | SUBROUTINE tide_update ( kt, idx, dta, td, jit, time_offset ) |
---|
270 | !!---------------------------------------------------------------------- |
---|
271 | !! *** SUBROUTINE tide_update *** |
---|
272 | !! |
---|
273 | !! ** Purpose : - Add tidal forcing to ssh, u2d and v2d OBC data arrays. |
---|
274 | !! |
---|
275 | !!---------------------------------------------------------------------- |
---|
276 | INTEGER, INTENT( in ) :: kt ! Main timestep counter |
---|
277 | !!gm doctor jit ==> kit |
---|
278 | TYPE(OBC_INDEX), INTENT( in ) :: idx ! OBC indices |
---|
279 | TYPE(OBC_DATA), INTENT(inout) :: dta ! OBC external data |
---|
280 | TYPE(TIDES_DATA),INTENT( in ) :: td ! tidal harmonics data |
---|
281 | INTEGER,INTENT(in),OPTIONAL :: jit ! Barotropic timestep counter (for timesplitting option) |
---|
282 | INTEGER,INTENT( in ), OPTIONAL :: time_offset ! time offset in units of timesteps. NB. if jit |
---|
283 | ! is present then units = subcycle timesteps. |
---|
284 | ! time_offset = 0 => get data at "now" time level |
---|
285 | ! time_offset = -1 => get data at "before" time level |
---|
286 | ! time_offset = +1 => get data at "after" time level |
---|
287 | ! etc. |
---|
288 | !! |
---|
289 | INTEGER :: itide, igrd, ib ! dummy loop indices |
---|
290 | INTEGER :: time_add ! time offset in units of timesteps |
---|
291 | REAL(wp) :: z_arg, z_sarg |
---|
292 | REAL(wp), DIMENSION(jptides_max) :: z_sist, z_cost |
---|
293 | !!---------------------------------------------------------------------- |
---|
294 | |
---|
295 | time_add = 0 |
---|
296 | IF( PRESENT(time_offset) ) THEN |
---|
297 | time_add = time_offset |
---|
298 | ENDIF |
---|
299 | |
---|
300 | ! Note tide phase speeds are in deg/hour, so we need to convert the |
---|
301 | ! elapsed time in seconds to hours by dividing by 3600.0 |
---|
302 | IF( PRESENT(jit) ) THEN |
---|
303 | z_arg = ( (kt-1) * rdt + (jit+time_add) * rdt / REAL(nn_baro,wp) ) * rad / 3600.0 |
---|
304 | ELSE |
---|
305 | z_arg = (kt+time_add) * rdt * rad / 3600.0 |
---|
306 | ENDIF |
---|
307 | |
---|
308 | DO itide = 1, td%ncpt |
---|
309 | z_sarg = z_arg * td%speed(itide) |
---|
310 | z_cost(itide) = COS( z_sarg ) |
---|
311 | z_sist(itide) = SIN( z_sarg ) |
---|
312 | END DO |
---|
313 | |
---|
314 | DO itide = 1, td%ncpt |
---|
315 | igrd=1 ! SSH on tracer grid. |
---|
316 | DO ib = 1, idx%nblenrim(igrd) |
---|
317 | dta%ssh(ib) = dta%ssh(ib) + td%ssh(ib,itide,1)*z_cost(itide) + td%ssh(ib,itide,2)*z_sist(itide) |
---|
318 | ! if(lwp) write(numout,*) 'z', ib, itide, dta%ssh(ib), td%ssh(ib,itide,1),td%ssh(ib,itide,2) |
---|
319 | END DO |
---|
320 | igrd=2 ! U grid |
---|
321 | DO ib=1, idx%nblenrim(igrd) |
---|
322 | dta%u2d(ib) = dta%u2d(ib) + td%u(ib,itide,1)*z_cost(itide) + td%u(ib,itide,2)*z_sist(itide) |
---|
323 | ! if(lwp) write(numout,*) 'u',ib,itide,utide(ib), td%u(ib,itide,1),td%u(ib,itide,2) |
---|
324 | END DO |
---|
325 | igrd=3 ! V grid |
---|
326 | DO ib=1, idx%nblenrim(igrd) |
---|
327 | dta%v2d(ib) = dta%v2d(ib) + td%v(ib,itide,1)*z_cost(itide) + td%v(ib,itide,2)*z_sist(itide) |
---|
328 | ! if(lwp) write(numout,*) 'v',ib,itide,vtide(ib), td%v(ib,itide,1),td%v(ib,itide,2) |
---|
329 | END DO |
---|
330 | END DO |
---|
331 | ! |
---|
332 | END SUBROUTINE tide_update |
---|
333 | |
---|
334 | !!gm doctor naming of dummy argument variables!!! and all local variables |
---|
335 | |
---|
336 | SUBROUTINE uvset( ihs, iday, imnth, iyr, f, z_vplu ) |
---|
337 | !!---------------------------------------------------------------------- |
---|
338 | !! *** SUBROUTINE uvset *** |
---|
339 | !! |
---|
340 | !! ** Purpose : - adjust tidal forcing for date factors |
---|
341 | !! |
---|
342 | !!---------------------------------------------------------------------- |
---|
343 | implicit none |
---|
344 | INTEGER, INTENT( in ) :: ihs ! Start time hours |
---|
345 | INTEGER, INTENT( in ) :: iday ! start time days |
---|
346 | INTEGER, INTENT( in ) :: imnth ! start time month |
---|
347 | INTEGER, INTENT( in ) :: iyr ! start time year |
---|
348 | !! |
---|
349 | !!gm nc is jptides_max ???? |
---|
350 | INTEGER , PARAMETER :: nc =15 ! maximum number of constituents |
---|
351 | CHARACTER(len=8), DIMENSION(nc) :: cname |
---|
352 | INTEGER :: year, vd, ivdy, ndc, i, k |
---|
353 | REAL(wp) :: ss, h, p, en, p1, rtd |
---|
354 | REAL(wp), DIMENSION(nc) :: f ! nodal correction |
---|
355 | REAL(wp), DIMENSION(nc) :: z_vplu ! phase correction |
---|
356 | REAL(wp), DIMENSION(nc) :: u, v, zig |
---|
357 | !! |
---|
358 | DATA cname/ 'q1' , 'o1' , 'p1' , 's1' , 'k1' , & |
---|
359 | & '2n2' , 'mu2' , 'n2' , 'nu2' , 'm2' , & |
---|
360 | & 'l2' , 't2' , 's2' , 'k2' , 'm4' / |
---|
361 | DATA zig/ .2338507481, .2433518789, .2610826055, .2617993878, .2625161701, & |
---|
362 | & .4868657873, .4881373225, .4963669182, .4976384533, .5058680490, & |
---|
363 | & .5153691799, .5228820265, .5235987756, .5250323419, 1.011736098 / |
---|
364 | !!---------------------------------------------------------------------- |
---|
365 | ! |
---|
366 | ! ihs - start time gmt on ... |
---|
367 | ! iday/imnth/iyr - date e.g. 12/10/87 |
---|
368 | ! |
---|
369 | CALL vday(iday,imnth,iyr,ivdy) |
---|
370 | vd = ivdy |
---|
371 | ! |
---|
372 | !rp note change of year number for d. blackman shpen |
---|
373 | !rp if(iyr.ge.1000) year=iyr-1900 |
---|
374 | !rp if(iyr.lt.1000) year=iyr |
---|
375 | year = iyr |
---|
376 | ! |
---|
377 | !.....year = year of required data |
---|
378 | !.....vd = day of required data..set up for 0000gmt day year |
---|
379 | ndc = nc |
---|
380 | !.....ndc = number of constituents allowed |
---|
381 | !!gm use rpi ? |
---|
382 | rtd = 360.0 / 6.2831852 |
---|
383 | DO i = 1, ndc |
---|
384 | zig(i) = zig(i)*rtd |
---|
385 | ! sigo(i)= zig(i) |
---|
386 | END DO |
---|
387 | |
---|
388 | !!gm try to avoid RETURN in F90 |
---|
389 | IF( year == 0 ) RETURN |
---|
390 | |
---|
391 | CALL shpen( year, vd, ss, h , p , en, p1 ) |
---|
392 | CALL ufset( p , en, u , f ) |
---|
393 | CALL vset ( ss , h , p , en, p1, v ) |
---|
394 | ! |
---|
395 | DO k = 1, nc |
---|
396 | z_vplu(k) = v(k) + u(k) |
---|
397 | z_vplu(k) = z_vplu(k) + dble(ihs) * zig(k) |
---|
398 | DO WHILE( z_vplu(k) < 0 ) |
---|
399 | z_vplu(k) = z_vplu(k) + 360.0 |
---|
400 | END DO |
---|
401 | DO WHILE( z_vplu(k) > 360. ) |
---|
402 | z_vplu(k) = z_vplu(k) - 360.0 |
---|
403 | END DO |
---|
404 | END DO |
---|
405 | ! |
---|
406 | END SUBROUTINE uvset |
---|
407 | |
---|
408 | |
---|
409 | SUBROUTINE vday( iday, imnth, iy, ivdy ) |
---|
410 | !!---------------------------------------------------------------------- |
---|
411 | !! *** SUBROUTINE vday *** |
---|
412 | !! |
---|
413 | !! ** Purpose : - adjust tidal forcing for date factors |
---|
414 | !! |
---|
415 | !!---------------------------------------------------------------------- |
---|
416 | INTEGER, INTENT(in ) :: iday, imnth, iy ! ???? |
---|
417 | INTEGER, INTENT( out) :: ivdy ! ??? |
---|
418 | !! |
---|
419 | INTEGER :: iyr |
---|
420 | !!------------------------------------------------------------------------------ |
---|
421 | |
---|
422 | !!gm nday_year in day mode is the variable compiuted here, no? |
---|
423 | !!gm nday_year , & !: curent day counted from jan 1st of the current year |
---|
424 | |
---|
425 | !calculate day number in year from day/month/year |
---|
426 | if(imnth.eq.1) ivdy=iday |
---|
427 | if(imnth.eq.2) ivdy=iday+31 |
---|
428 | if(imnth.eq.3) ivdy=iday+59 |
---|
429 | if(imnth.eq.4) ivdy=iday+90 |
---|
430 | if(imnth.eq.5) ivdy=iday+120 |
---|
431 | if(imnth.eq.6) ivdy=iday+151 |
---|
432 | if(imnth.eq.7) ivdy=iday+181 |
---|
433 | if(imnth.eq.8) ivdy=iday+212 |
---|
434 | if(imnth.eq.9) ivdy=iday+243 |
---|
435 | if(imnth.eq.10) ivdy=iday+273 |
---|
436 | if(imnth.eq.11) ivdy=iday+304 |
---|
437 | if(imnth.eq.12) ivdy=iday+334 |
---|
438 | iyr=iy |
---|
439 | if(mod(iyr,4).eq.0.and.imnth.gt.2) ivdy=ivdy+1 |
---|
440 | if(mod(iyr,100).eq.0.and.imnth.gt.2) ivdy=ivdy-1 |
---|
441 | if(mod(iyr,400).eq.0.and.imnth.gt.2) ivdy=ivdy+1 |
---|
442 | ! |
---|
443 | END SUBROUTINE vday |
---|
444 | |
---|
445 | !!doctor norme for dummy arguments |
---|
446 | |
---|
447 | SUBROUTINE shpen( year, vd, s, h, p, en, p1 ) |
---|
448 | !!---------------------------------------------------------------------- |
---|
449 | !! *** SUBROUTINE shpen *** |
---|
450 | !! |
---|
451 | !! ** Purpose : - calculate astronomical arguments for tides |
---|
452 | !! this version from d. blackman 30 nove 1990 |
---|
453 | !! |
---|
454 | !!---------------------------------------------------------------------- |
---|
455 | !!gm add INTENT in, out or inout.... DOCTOR name.... |
---|
456 | !!gm please do not use variable name with a single letter: impossible to search in a code |
---|
457 | INTEGER :: year,vd |
---|
458 | REAL(wp) :: s,h,p,en,p1 |
---|
459 | !! |
---|
460 | INTEGER :: yr,ilc,icent,it,iday,ild,ipos,nn,iyd |
---|
461 | REAL(wp) :: cycle,t,td,delt(84),delta,deltat |
---|
462 | !! |
---|
463 | DATA delt /-5.04, -3.90, -2.87, -0.58, 0.71, 1.80, & |
---|
464 | & 3.08, 4.63, 5.86, 7.21, 8.58, 10.50, 12.10, & |
---|
465 | & 12.49, 14.41, 15.59, 15.81, 17.52, 19.01, 18.39, & |
---|
466 | & 19.55, 20.36, 21.01, 21.81, 21.76, 22.35, 22.68, & |
---|
467 | & 22.94, 22.93, 22.69, 22.94, 23.20, 23.31, 23.63, & |
---|
468 | & 23.47, 23.68, 23.62, 23.53, 23.59, 23.99, 23.80, & |
---|
469 | & 24.20, 24.99, 24.97, 25.72, 26.21, 26.37, 26.89, & |
---|
470 | & 27.68, 28.13, 28.94, 29.42, 29.66, 30.29, 30.96, & |
---|
471 | & 31.09, 31.59, 31.52, 31.92, 32.45, 32.91, 33.39, & |
---|
472 | & 33.80, 34.23, 34.73, 35.40, 36.14, 36.99, 37.87, & |
---|
473 | & 38.75, 39.70, 40.70, 41.68, 42.82, 43.96, 45.00, & |
---|
474 | & 45.98, 47.00, 48.03, 49.10, 50.10, 50.97, 51.81, & |
---|
475 | & 52.57 / |
---|
476 | !!---------------------------------------------------------------------- |
---|
477 | |
---|
478 | cycle = 360.0 |
---|
479 | ilc = 0 |
---|
480 | icent = year / 100 |
---|
481 | yr = year - icent * 100 |
---|
482 | t = icent - 20 |
---|
483 | ! |
---|
484 | ! for the following equations |
---|
485 | ! time origin is fixed at 00 hr of jan 1st,2000. |
---|
486 | ! see notes by cartwright |
---|
487 | ! |
---|
488 | !!gm old coding style, use CASE instead and avoid GOTO (obsolescence in fortran 90) |
---|
489 | !!gm obsol( 1): Arithmetic IF statement is used ===> remove this in Fortran 90 |
---|
490 | it = icent - 20 |
---|
491 | if (it) 1,2,2 |
---|
492 | 1 iday = it/4 -it |
---|
493 | go to 3 |
---|
494 | 2 iday = (it+3)/4 - it |
---|
495 | ! |
---|
496 | ! t is in julian century |
---|
497 | ! correction in gegorian calander where only century year divisible |
---|
498 | ! by 4 is leap year. |
---|
499 | ! |
---|
500 | 3 continue |
---|
501 | ! |
---|
502 | td = 0.0 |
---|
503 | ! |
---|
504 | !!gm obsol( 1): Arithmetic IF statement is used ===> remove this in Fortran 90 |
---|
505 | if (yr) 4,5,4 |
---|
506 | ! |
---|
507 | 4 iyd = 365*yr |
---|
508 | ild = (yr-1)/4 |
---|
509 | if((icent - (icent/4)*4) .eq. 0) ilc = 1 |
---|
510 | td = iyd + ild + ilc |
---|
511 | ! |
---|
512 | 5 td = td + iday + vd -1.0 - 0.5 |
---|
513 | t = t + (td/36525.0) |
---|
514 | ! |
---|
515 | ipos=year-1899 |
---|
516 | if (ipos .lt. 0) go to 7 |
---|
517 | if (ipos .gt. 83) go to 6 |
---|
518 | ! |
---|
519 | delta = (delt(ipos+1)+delt(ipos))/2.0 |
---|
520 | go to 7 |
---|
521 | ! |
---|
522 | 6 delta= (65.0-50.5)/20.0*(year-1980)+50.5 |
---|
523 | ! |
---|
524 | 7 deltat = delta * 1.0e-6 |
---|
525 | ! |
---|
526 | !!gm precision of the computation : example for s it should be replace by: |
---|
527 | !!gm s = 218.3165 + (481267.8813 - 0.0016*t)*t + 152.0*deltat ==> more precise modify the last digits results |
---|
528 | s = 218.3165 + 481267.8813*t - 0.0016*t*t + 152.0*deltat |
---|
529 | h = 280.4661 + 36000.7698 *t + 0.0003*t*t + 11.0*deltat |
---|
530 | p = 83.3535 + 4069.0139 *t - 0.0103*t*t + deltat |
---|
531 | en = 234.9555 + 1934.1363 *t - 0.0021*t*t + deltat |
---|
532 | p1 = 282.9384 + 1.7195 *t + 0.0005*t*t |
---|
533 | ! |
---|
534 | nn = s / cycle |
---|
535 | s = s - nn * cycle |
---|
536 | IF( s < 0.e0 ) s = s + cycle |
---|
537 | ! |
---|
538 | nn = h / cycle |
---|
539 | h = h - cycle * nn |
---|
540 | IF( h < 0.e0 ) h = h + cycle |
---|
541 | ! |
---|
542 | nn = p / cycle |
---|
543 | p = p - cycle * nn |
---|
544 | IF( p < 0.e0) p = p + cycle |
---|
545 | ! |
---|
546 | nn = en / cycle |
---|
547 | en = en - cycle * nn |
---|
548 | IF( en < 0.e0 ) en = en + cycle |
---|
549 | en = cycle - en |
---|
550 | ! |
---|
551 | nn = p1 / cycle |
---|
552 | p1 = p1 - nn * cycle |
---|
553 | ! |
---|
554 | END SUBROUTINE shpen |
---|
555 | |
---|
556 | |
---|
557 | SUBROUTINE ufset( p, cn, b, a ) |
---|
558 | !!---------------------------------------------------------------------- |
---|
559 | !! *** SUBROUTINE ufset *** |
---|
560 | !! |
---|
561 | !! ** Purpose : - calculate nodal parameters for the tides |
---|
562 | !! |
---|
563 | !!---------------------------------------------------------------------- |
---|
564 | !!gm doctor naming of dummy argument variables!!! and all local variables |
---|
565 | !!gm nc is jptides_max ???? |
---|
566 | integer nc |
---|
567 | parameter (nc=15) |
---|
568 | REAL(wp) p,cn |
---|
569 | !! |
---|
570 | |
---|
571 | !!gm rad is already a public variable defined in phycst.F90 .... ==> doctor norme local real start with "z" |
---|
572 | REAL(wp) :: w1, w2, w3, w4, w5, w6, w7, w8, nw, pw, rad |
---|
573 | REAL(wp) :: a(nc), b(nc) |
---|
574 | INTEGER :: ndc, k |
---|
575 | !!---------------------------------------------------------------------- |
---|
576 | |
---|
577 | ndc = nc |
---|
578 | |
---|
579 | ! a=f , b =u |
---|
580 | ! t is zero as compared to tifa. |
---|
581 | !! use rad defined in phycst (i.e. add a USE phycst at the begining of the module |
---|
582 | rad = 6.2831852d0/360.0 |
---|
583 | pw = p * rad |
---|
584 | nw = cn * rad |
---|
585 | w1 = cos( nw ) |
---|
586 | w2 = cos( 2*nw ) |
---|
587 | w3 = cos( 3*nw ) |
---|
588 | w4 = sin( nw ) |
---|
589 | w5 = sin( 2*nw ) |
---|
590 | w6 = sin( 3*nw ) |
---|
591 | w7 = 1. - 0.2505 * COS( 2*pw ) - 0.1102 * COS(2*pw-nw ) & |
---|
592 | & - 0.156 * COS( 2*pw-2*nw ) - 0.037 * COS( nw ) |
---|
593 | w8 = - 0.2505 * SIN( 2*pw ) - 0.1102 * SIN(2*pw-nw ) & |
---|
594 | & - 0.0156 * SIN( 2*pw-2*nw ) - 0.037 * SIN( nw ) |
---|
595 | ! |
---|
596 | a(1) = 1.0089 + 0.1871 * w1 - 0.0147 * w2 + 0.0014 * w3 |
---|
597 | b(1) = 0.1885 * w4 - 0.0234 * w5 + 0.0033 * w6 |
---|
598 | ! q1 |
---|
599 | a(2) = a(1) |
---|
600 | b(2) = b(1) |
---|
601 | ! o1 |
---|
602 | a(3) = 1.0 |
---|
603 | b(3) = 0.0 |
---|
604 | ! p1 |
---|
605 | a(4) = 1.0 |
---|
606 | b(4) = 0.0 |
---|
607 | ! s1 |
---|
608 | a(5) = 1.0060+0.1150*w1- 0.0088*w2 +0.0006*w3 |
---|
609 | b(5) = -0.1546*w4 + 0.0119*w5 -0.0012*w6 |
---|
610 | ! k1 |
---|
611 | a(6) =1.0004 -0.0373*w1+ 0.0002*w2 |
---|
612 | b(6) = -0.0374*w4 |
---|
613 | ! 2n2 |
---|
614 | a(7) = a(6) |
---|
615 | b(7) = b(6) |
---|
616 | ! mu2 |
---|
617 | a(8) = a(6) |
---|
618 | b(8) = b(6) |
---|
619 | ! n2 |
---|
620 | a(9) = a(6) |
---|
621 | b(9) = b(6) |
---|
622 | ! nu2 |
---|
623 | a(10) = a(6) |
---|
624 | b(10) = b(6) |
---|
625 | ! m2 |
---|
626 | a(11) = SQRT( w7 * w7 + w8 * w8 ) |
---|
627 | b(11) = ATAN( w8 / w7 ) |
---|
628 | !!gmuse rpi instead of 3.141992 ??? true pi is rpi=3.141592653589793_wp ..... ???? |
---|
629 | IF( w7 < 0.e0 ) b(11) = b(11) + 3.141992 |
---|
630 | ! l2 |
---|
631 | a(12) = 1.0 |
---|
632 | b(12) = 0.0 |
---|
633 | ! t2 |
---|
634 | a(13)= a(12) |
---|
635 | b(13)= b(12) |
---|
636 | ! s2 |
---|
637 | a(14) = 1.0241+0.2863*w1+0.0083*w2 -0.0015*w3 |
---|
638 | b(14) = -0.3096*w4 + 0.0119*w5 - 0.0007*w6 |
---|
639 | ! k2 |
---|
640 | a(15) = a(6)*a(6) |
---|
641 | b(15) = 2*b(6) |
---|
642 | ! m4 |
---|
643 | !!gm old coding, remove GOTO and label of lines |
---|
644 | !!gm obsol( 1): Arithmetic IF statement is used ===> remove this in Fortran 90 |
---|
645 | DO 40 k = 1,ndc |
---|
646 | b(k) = b(k)/rad |
---|
647 | 32 if (b(k)) 34,35,35 |
---|
648 | 34 b(k) = b(k) + 360.0 |
---|
649 | go to 32 |
---|
650 | 35 if (b(k)-360.0) 40,37,37 |
---|
651 | 37 b(k) = b(k)-360.0 |
---|
652 | go to 35 |
---|
653 | 40 continue |
---|
654 | ! |
---|
655 | END SUBROUTINE ufset |
---|
656 | |
---|
657 | |
---|
658 | SUBROUTINE vset( s,h,p,en,p1,v) |
---|
659 | !!---------------------------------------------------------------------- |
---|
660 | !! *** SUBROUTINE vset *** |
---|
661 | !! |
---|
662 | !! ** Purpose : - calculate tidal phases for 0000gmt on start day of run |
---|
663 | !! |
---|
664 | !!---------------------------------------------------------------------- |
---|
665 | !!gm doctor naming of dummy argument variables!!! and all local variables |
---|
666 | !!gm nc is jptides_max ???? |
---|
667 | !!gm en argument is not used: suppress it ? |
---|
668 | integer nc |
---|
669 | parameter (nc=15) |
---|
670 | real(wp) s,h,p,en,p1 |
---|
671 | real(wp) v(nc) |
---|
672 | !! |
---|
673 | integer ndc, k |
---|
674 | !!---------------------------------------------------------------------- |
---|
675 | |
---|
676 | ndc = nc |
---|
677 | ! v s are computed here. |
---|
678 | v(1) =-3*s +h +p +270 ! Q1 |
---|
679 | v(2) =-2*s +h +270.0 ! O1 |
---|
680 | v(3) =-h +270 ! P1 |
---|
681 | v(4) =180 ! S1 |
---|
682 | v(5) =h +90.0 ! K1 |
---|
683 | v(6) =-4*s +2*h +2*p ! 2N2 |
---|
684 | v(7) =-4*(s-h) ! MU2 |
---|
685 | v(8) =-3*s +2*h +p ! N2 |
---|
686 | v(9) =-3*s +4*h -p ! MU2 |
---|
687 | v(10) =-2*s +2*h ! M2 |
---|
688 | v(11) =-s +2*h -p +180 ! L2 |
---|
689 | v(12) =-h +p1 ! T2 |
---|
690 | v(13) =0 ! S2 |
---|
691 | v(14) =h+h ! K2 |
---|
692 | v(15) =2*v(10) ! M4 |
---|
693 | ! |
---|
694 | !!gm old coding, remove GOTO and label of lines |
---|
695 | !!gm obsol( 1): Arithmetic IF statement is used ===> remove this in Fortran 90 |
---|
696 | do 72 k = 1, ndc |
---|
697 | 69 if( v(k) ) 70,71,71 |
---|
698 | 70 v(k) = v(k)+360.0 |
---|
699 | go to 69 |
---|
700 | 71 if( v(k) - 360.0 ) 72,73,73 |
---|
701 | 73 v(k) = v(k)-360.0 |
---|
702 | go to 71 |
---|
703 | 72 continue |
---|
704 | ! |
---|
705 | END SUBROUTINE vset |
---|
706 | |
---|
707 | #else |
---|
708 | !!---------------------------------------------------------------------- |
---|
709 | !! Dummy module NO Unstruct Open Boundary Conditions for tides |
---|
710 | !!---------------------------------------------------------------------- |
---|
711 | !!gm are you sure we need to define filtide and tide_cpt ? |
---|
712 | CHARACTER(len=80), PUBLIC :: filtide !: Filename root for tidal input files |
---|
713 | CHARACTER(len=4 ), PUBLIC, DIMENSION(1) :: tide_cpt !: Names of tidal components used. |
---|
714 | |
---|
715 | CONTAINS |
---|
716 | SUBROUTINE tide_init ! Empty routine |
---|
717 | END SUBROUTINE tide_init |
---|
718 | SUBROUTINE tide_data ! Empty routine |
---|
719 | END SUBROUTINE tide_data |
---|
720 | SUBROUTINE tide_update( kt, kit ) ! Empty routine |
---|
721 | WRITE(*,*) 'tide_update: You should not have seen this print! error?', kt, kit |
---|
722 | END SUBROUTINE tide_update |
---|
723 | #endif |
---|
724 | |
---|
725 | !!====================================================================== |
---|
726 | END MODULE bdytides |
---|