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 | !! 3.4 ! 2011 (D. Storkey) rewrite in preparation for OBC-BDY merge |
---|
11 | !! 3.4 ! 2013 (J. Harle) rewite to used tide_mod for phase and nodal |
---|
12 | !! corrections every day |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | #if defined key_bdy |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! 'key_bdy' Open Boundary Condition |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | !! PUBLIC |
---|
19 | !! tide_init : read of namelist and initialisation of tidal harmonics data |
---|
20 | !! tide_update : calculation of tidal forcing at each timestep |
---|
21 | !!---------------------------------------------------------------------- |
---|
22 | USE timing ! Timing |
---|
23 | USE oce ! ocean dynamics and tracers |
---|
24 | USE dom_oce ! ocean space and time domain |
---|
25 | USE iom |
---|
26 | USE in_out_manager ! I/O units |
---|
27 | USE phycst ! physical constants |
---|
28 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
29 | USE bdy_par ! Unstructured boundary parameters |
---|
30 | USE bdy_oce ! ocean open boundary conditions |
---|
31 | USE fldread, ONLY: fld_map |
---|
32 | USE daymod ! calendar |
---|
33 | USE tide_mod |
---|
34 | USE ioipsl, ONLY : ymds2ju ! for calendar |
---|
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 | REAL(wp), POINTER, DIMENSION(:,:,:) :: sshr !: Tidal constituents : SSH (reference) |
---|
49 | REAL(wp), POINTER, DIMENSION(:,:,:) :: ur !: Tidal constituents : U (reference) |
---|
50 | REAL(wp), POINTER, DIMENSION(:,:,:) :: vr !: Tidal constituents : V (reference) |
---|
51 | END TYPE TIDES_DATA |
---|
52 | |
---|
53 | TYPE(TIDES_DATA), PUBLIC, DIMENSION(jp_bdy), TARGET :: tides !: External tidal harmonics data |
---|
54 | |
---|
55 | INTEGER, ALLOCATABLE, DIMENSION(:) :: bdy_ntide |
---|
56 | REAL(wp), ALLOCATABLE, DIMENSION(:) :: bdy_omega_tide |
---|
57 | REAL(wp), ALLOCATABLE, DIMENSION(:) :: bdy_v0tide, & |
---|
58 | bdy_blank, & |
---|
59 | bdy_utide, & |
---|
60 | bdy_ftide, & |
---|
61 | rbdy_ftide |
---|
62 | LOGICAL :: ln_tide_date !: =T correct tide phases and amplitude for model start date |
---|
63 | LOGICAL :: ln_tide_v0 !: =T correct tide phases and amplitude for model start date |
---|
64 | INTEGER :: nn_tide_date !: yyyymmdd reference date of tidal data |
---|
65 | INTEGER :: bdy_nn_tide |
---|
66 | INTEGER :: bdy_kt_tide ! Main tide timestep counter |
---|
67 | INTEGER :: bdy_tide_offset ! Main tide timestep counter |
---|
68 | |
---|
69 | !!---------------------------------------------------------------------- |
---|
70 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
71 | !! $Id$ |
---|
72 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
73 | !!---------------------------------------------------------------------- |
---|
74 | CONTAINS |
---|
75 | |
---|
76 | SUBROUTINE tide_init |
---|
77 | !!---------------------------------------------------------------------- |
---|
78 | !! *** SUBROUTINE tide_init *** |
---|
79 | !! |
---|
80 | !! ** Purpose : - Read in namelist for tides and initialise external |
---|
81 | !! tidal harmonics data |
---|
82 | !! |
---|
83 | !!---------------------------------------------------------------------- |
---|
84 | !! namelist variables |
---|
85 | !!------------------- |
---|
86 | CHARACTER(len=80) :: filtide !: Filename root for tidal input files |
---|
87 | CHARACTER(len= 4), DIMENSION(jpmax_harmo) :: tide_cpt !: Names of tidal components used. |
---|
88 | !! |
---|
89 | INTEGER :: ib_bdy, itide, ib, ji !: dummy loop indices |
---|
90 | INTEGER :: inum, igrd |
---|
91 | INTEGER :: lcl_ryear, lcl_rmonth, lcl_rday |
---|
92 | INTEGER, DIMENSION(3) :: ilen0 !: length of boundary data (from OBC arrays) |
---|
93 | CHARACTER(len=80) :: clfile !: full file name for tidal input file |
---|
94 | REAL(wp) :: z_arg, z_atde, z_btde, z1t, z2t, fdayn, fdayr |
---|
95 | REAL(wp),ALLOCATABLE, DIMENSION(:,:,:) :: dta_read !: work space to read in tidal harmonics data |
---|
96 | !! |
---|
97 | TYPE(TIDES_DATA), POINTER :: td !: local short cut |
---|
98 | !! |
---|
99 | NAMELIST/nambdy_tide/filtide, tide_cpt, ln_tide_date, nn_tide_date, ln_tide_v0 |
---|
100 | !!---------------------------------------------------------------------- |
---|
101 | |
---|
102 | IF( nn_timing == 1 ) CALL timing_start('tide_init') |
---|
103 | |
---|
104 | IF(lwp) WRITE(numout,*) |
---|
105 | IF(lwp) WRITE(numout,*) 'tide_init : initialization of tidal harmonic forcing at open boundaries' |
---|
106 | IF(lwp) WRITE(numout,*) '~~~~~~~~~' |
---|
107 | |
---|
108 | REWIND(numnam) |
---|
109 | DO ib_bdy = 1, nb_bdy |
---|
110 | IF( nn_dyn2d_dta(ib_bdy) .ge. 2 ) THEN |
---|
111 | |
---|
112 | td => tides(ib_bdy) |
---|
113 | |
---|
114 | ! Namelist nambdy_tide : tidal harmonic forcing at open boundaries |
---|
115 | ln_tide_date = .false. |
---|
116 | ln_tide_v0 = .false. |
---|
117 | filtide(:) = '' |
---|
118 | tide_cpt(:) = '' |
---|
119 | |
---|
120 | ! Initialise bdy_ky_tide: updated in tide_update if using time correction otherwise defaults to 1 |
---|
121 | bdy_kt_tide=1 |
---|
122 | |
---|
123 | ! Don't REWIND here - may need to read more than one of these namelists. |
---|
124 | READ ( numnam, nambdy_tide ) |
---|
125 | ! ! Count number of components specified |
---|
126 | td%ncpt = 0 |
---|
127 | DO itide = 1, jpmax_harmo |
---|
128 | IF( tide_cpt(itide) /= '' ) THEN |
---|
129 | td%ncpt = td%ncpt + 1 |
---|
130 | ENDIF |
---|
131 | END DO |
---|
132 | |
---|
133 | CALL tide_init_Wave |
---|
134 | |
---|
135 | ! Find constituents in standard list |
---|
136 | ALLOCATE(bdy_ntide (td%ncpt)) |
---|
137 | |
---|
138 | DO itide=1,td%ncpt |
---|
139 | bdy_ntide(itide)=0 |
---|
140 | DO ji=1,jpmax_harmo |
---|
141 | IF ( TRIM( tide_cpt(itide) ) .eq. Wave(ji)%cname_tide) THEN |
---|
142 | bdy_ntide(itide) = ji |
---|
143 | EXIT |
---|
144 | END IF |
---|
145 | END DO |
---|
146 | IF (bdy_ntide(itide).eq.0) THEN |
---|
147 | CALL ctl_stop( 'BDYTIDE tidal components do not match up with tide.h90' ) |
---|
148 | ENDIF |
---|
149 | END DO |
---|
150 | |
---|
151 | ! Fill in phase speeds from tide_pulse |
---|
152 | ALLOCATE(bdy_omega_tide(td%ncpt)) |
---|
153 | CALL tide_pulse( bdy_omega_tide, bdy_ntide ,td%ncpt) |
---|
154 | |
---|
155 | ALLOCATE( td%speed(td%ncpt) ) |
---|
156 | td%speed = bdy_omega_tide(1:td%ncpt) |
---|
157 | |
---|
158 | ! ! Parameter control and print |
---|
159 | IF( td%ncpt < 1 ) THEN |
---|
160 | CALL ctl_stop( ' Did not find any tidal components in namelist nambdy_tide' ) |
---|
161 | ELSE |
---|
162 | IF(lwp) WRITE(numout,*) ' Namelist nambdy_tide : tidal harmonic forcing at open boundaries' |
---|
163 | IF(lwp) WRITE(numout,*) ' tidal components specified ', td%ncpt |
---|
164 | IF(lwp) WRITE(numout,*) ' ', tide_cpt(1:td%ncpt) |
---|
165 | IF(lwp) WRITE(numout,*) ' associated phase speeds (deg/hr) : ' |
---|
166 | IF(lwp) WRITE(numout,*) ' ', td%speed(1:td%ncpt) |
---|
167 | ENDIF |
---|
168 | |
---|
169 | ! Allocate space for tidal harmonics data - |
---|
170 | ! get size from OBC data arrays |
---|
171 | ! --------------------------------------- |
---|
172 | |
---|
173 | ilen0(1) = SIZE( dta_bdy(ib_bdy)%ssh ) |
---|
174 | ALLOCATE( td%ssh( ilen0(1), td%ncpt, 2 ) ) |
---|
175 | ALLOCATE( td%sshr( ilen0(1), td%ncpt, 2 ) ) |
---|
176 | |
---|
177 | ilen0(2) = SIZE( dta_bdy(ib_bdy)%u2d ) |
---|
178 | ALLOCATE( td%u( ilen0(2), td%ncpt, 2 ) ) |
---|
179 | ALLOCATE( td%ur( ilen0(2), td%ncpt, 2 ) ) |
---|
180 | |
---|
181 | ilen0(3) = SIZE( dta_bdy(ib_bdy)%v2d ) |
---|
182 | ALLOCATE( td%v( ilen0(3), td%ncpt, 2 ) ) |
---|
183 | ALLOCATE( td%vr( ilen0(3), td%ncpt, 2 ) ) |
---|
184 | |
---|
185 | ALLOCATE( dta_read( MAXVAL(ilen0), 1, 1 ) ) |
---|
186 | |
---|
187 | ! Set day length in timesteps for use if making phase and nodal corrections |
---|
188 | bdy_nn_tide=NINT(rday/rdt) |
---|
189 | |
---|
190 | |
---|
191 | ALLOCATE(bdy_v0tide (td%ncpt)) |
---|
192 | ALLOCATE(bdy_blank (td%ncpt)) |
---|
193 | ALLOCATE(bdy_utide (td%ncpt)) |
---|
194 | ALLOCATE(bdy_ftide (td%ncpt)) |
---|
195 | ALLOCATE(rbdy_ftide (td%ncpt)) |
---|
196 | |
---|
197 | ! Open files and read in tidal forcing data |
---|
198 | ! ----------------------------------------- |
---|
199 | |
---|
200 | DO itide = 1, td%ncpt |
---|
201 | ! ! SSH fields |
---|
202 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_T.nc' |
---|
203 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
204 | CALL iom_open( clfile, inum ) |
---|
205 | CALL fld_map( inum, 'z1' , dta_read(1:ilen0(1),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,1) ) |
---|
206 | td%ssh(:,itide,1) = dta_read(1:ilen0(1),1,1) |
---|
207 | CALL fld_map( inum, 'z2' , dta_read(1:ilen0(1),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,1) ) |
---|
208 | td%ssh(:,itide,2) = dta_read(1:ilen0(1),1,1) |
---|
209 | CALL iom_close( inum ) |
---|
210 | ! ! U fields |
---|
211 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_U.nc' |
---|
212 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
213 | CALL iom_open( clfile, inum ) |
---|
214 | CALL fld_map( inum, 'u1' , dta_read(1:ilen0(2),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,2) ) |
---|
215 | td%u(:,itide,1) = dta_read(1:ilen0(2),1,1) |
---|
216 | CALL fld_map( inum, 'u2' , dta_read(1:ilen0(2),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,2) ) |
---|
217 | td%u(:,itide,2) = dta_read(1:ilen0(2),1,1) |
---|
218 | CALL iom_close( inum ) |
---|
219 | ! ! V fields |
---|
220 | clfile = TRIM(filtide)//TRIM(tide_cpt(itide))//'_grid_V.nc' |
---|
221 | IF(lwp) WRITE(numout,*) 'Reading data from file ', clfile |
---|
222 | CALL iom_open( clfile, inum ) |
---|
223 | CALL fld_map( inum, 'v1' , dta_read(1:ilen0(3),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,3) ) |
---|
224 | td%v(:,itide,1) = dta_read(1:ilen0(3),1,1) |
---|
225 | CALL fld_map( inum, 'v2' , dta_read(1:ilen0(3),1:1,1:1) , 1, idx_bdy(ib_bdy)%nbmap(:,3) ) |
---|
226 | td%v(:,itide,2) = dta_read(1:ilen0(3),1,1) |
---|
227 | CALL iom_close( inum ) |
---|
228 | ! |
---|
229 | END DO ! end loop on tidal components |
---|
230 | |
---|
231 | IF( ln_tide_date .and. ln_tide_v0 ) THEN ! correct for date factors: gather v0 |
---|
232 | CALL tide_harmo(bdy_omega_tide, bdy_v0tide, bdy_utide, bdy_ftide, bdy_ntide, td%ncpt, nn_tide_date) |
---|
233 | |
---|
234 | lcl_ryear = INT(nn_tide_date / 10000 ) |
---|
235 | lcl_rmonth = INT((nn_tide_date - lcl_ryear * 10000 ) / 100 ) |
---|
236 | lcl_rday = INT(nn_tide_date - lcl_ryear * 10000 - lcl_rmonth * 100) |
---|
237 | nyear = int(ndate0 / 10000 ) ! initial year |
---|
238 | nmonth = int((ndate0 - nyear * 10000 ) / 100 ) ! initial month |
---|
239 | nday = int(ndate0 - nyear * 10000 - nmonth * 100) |
---|
240 | CALL ymds2ju( nyear, nmonth, nday, 0._wp, fdayn ) |
---|
241 | CALL ymds2ju( lcl_ryear, lcl_rmonth, lcl_rday, 0._wp, fdayr ) |
---|
242 | bdy_tide_offset = NINT( fdayn - fdayr ) * 86400 |
---|
243 | IF(lwp) WRITE(numout,*) ' BDYTIDE offset ' |
---|
244 | IF(lwp) WRITE(numout,*) ' ', lcl_ryear, lcl_rmonth, lcl_rday |
---|
245 | IF(lwp) WRITE(numout,*) ' ', nyear, nmonth, nday |
---|
246 | IF(lwp) WRITE(numout,*) ' ', fdayn, fdayr, bdy_tide_offset |
---|
247 | ELSE |
---|
248 | bdy_v0tide(:)=0 |
---|
249 | bdy_utide(:)=0 |
---|
250 | bdy_ftide(:)=1 |
---|
251 | bdy_tide_offset = 0 |
---|
252 | IF(lwp) WRITE(numout,*) ' BDYTIDE offset ', bdy_tide_offset |
---|
253 | ENDIF |
---|
254 | |
---|
255 | ! Pass tidal forcing data to reference arrays for date correction to tidal harmonics |
---|
256 | |
---|
257 | DO itide = 1, td%ncpt ! loop on tidal components |
---|
258 | ! ! elevation |
---|
259 | igrd = 1 |
---|
260 | DO ib = 1, ilen0(igrd) |
---|
261 | td%sshr(ib,itide,1) = td%ssh(ib,itide,1) |
---|
262 | td%sshr(ib,itide,2) = td%ssh(ib,itide,2) |
---|
263 | END DO |
---|
264 | ! ! u |
---|
265 | igrd = 2 |
---|
266 | DO ib = 1, ilen0(igrd) |
---|
267 | td%ur(ib,itide,1) = td%u(ib,itide,1) |
---|
268 | td%ur(ib,itide,2) = td%u(ib,itide,2) |
---|
269 | END DO |
---|
270 | ! ! v |
---|
271 | igrd = 3 |
---|
272 | DO ib = 1, ilen0(igrd) |
---|
273 | td%vr(ib,itide,1) = td%v(ib,itide,1) |
---|
274 | td%vr(ib,itide,2) = td%v(ib,itide,2) |
---|
275 | ENDDO |
---|
276 | ENDDO ! loop on tidal components |
---|
277 | |
---|
278 | IF(lwp) WRITE(numout,*) 'BDYTIDE: summary of mappings' |
---|
279 | DO itide = 1, td%ncpt ! loop on tidal components |
---|
280 | IF(lwp) WRITE(numout,'(2i3,x,a)') itide, bdy_ntide(itide), tide_cpt(itide) |
---|
281 | ENDDO |
---|
282 | |
---|
283 | ! |
---|
284 | ENDIF ! nn_dyn2d_dta(ib_bdy) .ge. 2 |
---|
285 | ! |
---|
286 | END DO ! loop on ib_bdy |
---|
287 | |
---|
288 | IF( nn_timing == 1 ) CALL timing_stop('tide_init') |
---|
289 | |
---|
290 | END SUBROUTINE tide_init |
---|
291 | |
---|
292 | |
---|
293 | SUBROUTINE tide_update ( kt, idx, dta, td, jit, time_offset ) |
---|
294 | !!---------------------------------------------------------------------- |
---|
295 | !! *** SUBROUTINE tide_update *** |
---|
296 | !! |
---|
297 | !! ** Purpose : - Add tidal forcing to ssh, u2d and v2d OBC data arrays. |
---|
298 | !! |
---|
299 | !!---------------------------------------------------------------------- |
---|
300 | INTEGER, INTENT( in ) :: kt ! Main timestep counter |
---|
301 | !!gm doctor jit ==> kit |
---|
302 | TYPE(OBC_INDEX), INTENT( in ) :: idx ! OBC indices |
---|
303 | TYPE(OBC_DATA), INTENT(inout) :: dta ! OBC external data |
---|
304 | TYPE(TIDES_DATA),INTENT(inout) :: td ! tidal harmonics data |
---|
305 | INTEGER,INTENT(in),OPTIONAL :: jit ! Barotropic timestep counter (for timesplitting option) |
---|
306 | INTEGER,INTENT( in ), OPTIONAL :: time_offset ! time offset in units of timesteps. NB. if jit |
---|
307 | ! is present then units = subcycle timesteps. |
---|
308 | ! time_offset = 0 => get data at "now" time level |
---|
309 | ! time_offset = -1 => get data at "before" time level |
---|
310 | ! time_offset = +1 => get data at "after" time level |
---|
311 | ! etc. |
---|
312 | !! |
---|
313 | INTEGER :: itide, igrd, ib ! dummy loop indices |
---|
314 | INTEGER :: time_add ! time offset in units of timesteps |
---|
315 | INTEGER :: sub_step ! dummy for jit (probably not required as |
---|
316 | ! timesplitting always used?) |
---|
317 | REAL(wp) :: z_arg, z_sarg |
---|
318 | REAL(wp), DIMENSION(jpmax_harmo) :: z_sist, z_cost |
---|
319 | REAL(wp) :: z_atde, z_btde |
---|
320 | REAL(wp) :: z1t, z2t |
---|
321 | !!---------------------------------------------------------------------- |
---|
322 | |
---|
323 | IF( nn_timing == 1 ) CALL timing_start('tide_update') |
---|
324 | |
---|
325 | time_add = 0 |
---|
326 | IF( PRESENT(time_offset) ) THEN |
---|
327 | time_add = time_offset |
---|
328 | ENDIF |
---|
329 | |
---|
330 | ! Phase corrections for the current day |
---|
331 | |
---|
332 | sub_step = 1 |
---|
333 | IF( PRESENT(jit) ) THEN |
---|
334 | sub_step = jit |
---|
335 | ENDIF |
---|
336 | |
---|
337 | IF( ln_tide_date ) THEN ! correct for date factors |
---|
338 | |
---|
339 | IF ( ( MOD( kt - 1, bdy_nn_tide ) == 0 ) .and. (sub_step==1) ) THEN |
---|
340 | IF ( ln_tide_v0 ) THEN |
---|
341 | bdy_kt_tide = 1 |
---|
342 | CALL tide_harmo(bdy_omega_tide, bdy_blank, bdy_utide, bdy_ftide, bdy_ntide, td%ncpt, ndastp) |
---|
343 | ELSE |
---|
344 | bdy_kt_tide = kt |
---|
345 | CALL tide_harmo(bdy_omega_tide, bdy_v0tide, bdy_utide, bdy_ftide, bdy_ntide, td%ncpt, ndastp) |
---|
346 | ENDIF |
---|
347 | |
---|
348 | DO itide = 1, td%ncpt ! loop on tidal components |
---|
349 | IF(lwp) WRITE(numout,*) 'BDYTIDE CORR:', itide, bdy_omega_tide(itide), bdy_v0tide(itide), & |
---|
350 | & bdy_utide(itide), bdy_ftide(itide) |
---|
351 | ENDDO |
---|
352 | |
---|
353 | ! Make adjustment for reference date in tidal harmonic data |
---|
354 | IF(lwp) WRITE(numout,*) 'BDYTIDE: nodal and phase correction at the start of day ', & |
---|
355 | & (kt-1)*rdt/rday + 1 |
---|
356 | |
---|
357 | DO itide = 1, td%ncpt ! loop on tidal components |
---|
358 | z_arg = bdy_utide(itide)+bdy_v0tide(itide) |
---|
359 | z_atde= bdy_ftide(itide)* cos(z_arg) |
---|
360 | z_btde= bdy_ftide(itide)* sin(z_arg) |
---|
361 | ! ! elevation |
---|
362 | igrd = 1 |
---|
363 | DO ib = 1, idx%nblenrim(igrd) |
---|
364 | z1t = z_atde * td%sshr(ib,itide,1) + z_btde * td%sshr(ib,itide,2) |
---|
365 | z2t = z_atde * td%sshr(ib,itide,2) - z_btde * td%sshr(ib,itide,1) |
---|
366 | td%ssh(ib,itide,1) = z1t |
---|
367 | td%ssh(ib,itide,2) = z2t |
---|
368 | END DO |
---|
369 | ! ! u |
---|
370 | igrd = 2 |
---|
371 | DO ib = 1, idx%nblenrim(igrd) |
---|
372 | z1t = z_atde * td%ur(ib,itide,1) + z_btde * td%ur(ib,itide,2) |
---|
373 | z2t = z_atde * td%ur(ib,itide,2) - z_btde * td%ur(ib,itide,1) |
---|
374 | td%u(ib,itide,1) = z1t |
---|
375 | td%u(ib,itide,2) = z2t |
---|
376 | END DO |
---|
377 | ! ! v |
---|
378 | igrd = 3 |
---|
379 | DO ib = 1, idx%nblenrim(igrd) |
---|
380 | z1t = z_atde * td%vr(ib,itide,1) + z_btde * td%vr(ib,itide,2) |
---|
381 | z2t = z_atde * td%vr(ib,itide,2) - z_btde * td%vr(ib,itide,1) |
---|
382 | td%v(ib,itide,1) = z1t |
---|
383 | td%v(ib,itide,2) = z2t |
---|
384 | ENDDO |
---|
385 | ENDDO ! loop on tidal components |
---|
386 | |
---|
387 | ENDIF |
---|
388 | |
---|
389 | ENDIF ! correct for date factors |
---|
390 | |
---|
391 | IF( PRESENT(jit) ) THEN |
---|
392 | IF( ln_tide_date ) THEN ! correct for date factors |
---|
393 | z_arg = ( (kt-bdy_kt_tide) * rdt + bdy_tide_offset + (jit+time_add) * rdt / REAL(nn_baro,wp) ) |
---|
394 | ELSE |
---|
395 | z_arg = ( (kt-1) * rdt + (jit+time_add) * rdt / REAL(nn_baro,wp) ) |
---|
396 | ENDIF |
---|
397 | ELSE |
---|
398 | IF(lwp) WRITE(numout,*) 'BDYTIDE: should I be in here?' |
---|
399 | IF( ln_tide_date ) THEN ! correct for date factors |
---|
400 | z_arg = (kt+time_add-bdy_kt_tide+1) * rdt + bdy_tide_offset |
---|
401 | ELSE |
---|
402 | z_arg = (kt+time_add) * rdt |
---|
403 | ENDIF |
---|
404 | ENDIF |
---|
405 | |
---|
406 | DO itide = 1, td%ncpt |
---|
407 | z_sarg = z_arg * td%speed(itide) |
---|
408 | z_cost(itide) = COS( z_sarg ) |
---|
409 | z_sist(itide) = SIN( z_sarg ) |
---|
410 | END DO |
---|
411 | |
---|
412 | DO itide = 1, td%ncpt |
---|
413 | igrd=1 ! SSH on tracer grid. |
---|
414 | DO ib = 1, idx%nblenrim(igrd) |
---|
415 | dta%ssh(ib) = dta%ssh(ib) + td%ssh(ib,itide,1)*z_cost(itide) + td%ssh(ib,itide,2)*z_sist(itide) |
---|
416 | IF ( (idx%nbmap(ib,igrd) == 100 .and. (itide==10)) .and. (sub_step==1) ) THEN |
---|
417 | write(numout,*) 'z', ib, idx%nbmap(ib,igrd), idx%nbi(ib,igrd), idx%nbj(ib,igrd), & |
---|
418 | & itide, (td%ssh(ib,itide,1)*z_cost(itide) + td%ssh(ib,itide,2)*z_sist(itide)) |
---|
419 | ENDIF |
---|
420 | ! IF ( (idx%nbmap(ib,igrd) == 100 .and. (itide==10)) ) THEN |
---|
421 | ! write(numout,*) 'z', ib, idx%nbmap(ib,igrd), idx%nbi(ib,igrd), idx%nbj(ib,igrd), & |
---|
422 | ! & itide, (td%ssh(ib,itide,1)*z_cost(itide) - td%ssh(ib,itide,2)*z_sist(itide)) |
---|
423 | ! ENDIF |
---|
424 | END DO |
---|
425 | igrd=2 ! U grid |
---|
426 | DO ib=1, idx%nblenrim(igrd) |
---|
427 | dta%u2d(ib) = dta%u2d(ib) + td%u(ib,itide,1)*z_cost(itide) + td%u(ib,itide,2)*z_sist(itide) |
---|
428 | ! if(lwp) write(numout,*) 'u',ib,itide,utide(ib), td%u(ib,itide,1),td%u(ib,itide,2) |
---|
429 | END DO |
---|
430 | igrd=3 ! V grid |
---|
431 | DO ib=1, idx%nblenrim(igrd) |
---|
432 | dta%v2d(ib) = dta%v2d(ib) + td%v(ib,itide,1)*z_cost(itide) + td%v(ib,itide,2)*z_sist(itide) |
---|
433 | ! if(lwp) write(numout,*) 'v',ib,itide,vtide(ib), td%v(ib,itide,1),td%v(ib,itide,2) |
---|
434 | END DO |
---|
435 | END DO |
---|
436 | ! |
---|
437 | IF( nn_timing == 1 ) CALL timing_stop('tide_update') |
---|
438 | ! |
---|
439 | END SUBROUTINE tide_update |
---|
440 | |
---|
441 | #else |
---|
442 | !!---------------------------------------------------------------------- |
---|
443 | !! Dummy module NO Unstruct Open Boundary Conditions for tides |
---|
444 | !!---------------------------------------------------------------------- |
---|
445 | !!gm are you sure we need to define filtide and tide_cpt ? |
---|
446 | CHARACTER(len=80), PUBLIC :: filtide !: Filename root for tidal input files |
---|
447 | CHARACTER(len=4 ), PUBLIC, DIMENSION(1) :: tide_cpt !: Names of tidal components used. |
---|
448 | |
---|
449 | CONTAINS |
---|
450 | SUBROUTINE tide_init ! Empty routine |
---|
451 | END SUBROUTINE tide_init |
---|
452 | SUBROUTINE tide_data ! Empty routine |
---|
453 | END SUBROUTINE tide_data |
---|
454 | SUBROUTINE tide_update( kt, kit ) ! Empty routine |
---|
455 | WRITE(*,*) 'tide_update: You should not have seen this print! error?', kt, kit |
---|
456 | END SUBROUTINE tide_update |
---|
457 | #endif |
---|
458 | |
---|
459 | !!====================================================================== |
---|
460 | END MODULE bdytides |
---|