ROMS
Loading...
Searching...
No Matches
mod_tides.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 MODULE mod_tides
3#if defined SSH_TIDES || defined UV_TIDES
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! Tidal Components: !
13! !
14! Each of the following arrays has a dimension in tidal components !
15! classified by period: !
16! !
17! semi-diurnal: M2, S2, N2, K2 (12.42, 12.00, 12.66, 11.97h) !
18! diurnal: K1, O1, P1, Q1 (23.93, 25.82, 24.07, 26.87h) !
19! !
20! and other longer periods. The order of these tidal components is !
21! irrelevant here. The number of components to USE is depends on !
22! the regional application. !
23! !
24! CosOmega Cosine tidal harmonics for current omega(t). !
25! SinOmega Sine tidal harmonics for current omega(t). !
26! SSH_Tamp Tidal elevation amplitude (m) at RHO-points. !
27! SSH_Tphase Tidal elevation phase (degrees/360) at RHO-points. !
28! Tperiod Tidal period (s). !
29! UV_Tangle Tidal current angle (radians; counterclockwise !
30! from EAST and rotated to curvilinear grid) at !
31! RHO-points. !
32! UV_Tmajor Maximum tidal current: tidal ellipse major axis !
33! (m/s) at RHO-points. !
34! UV_Tminor Minimum tidal current: tidal ellipse minor axis !
35! (m/s) at RHO-points. !
36! UV_Tphase Tidal current phase (degrees/360) at RHO-points. !
37! !
38# if defined AVERAGES && defined AVERAGES_DETIDE
39! !
40! Detided time-averaged fields via least-squares fitting. Notice that !
41! the harmonics for the state variable have an extra dimension of !
42! size (0:2*NTC) to store several terms: !
43! !
44! index 0 mean term (accumulated sum) !
45! 1:NTC accumulated sine terms !
46! NTC+1:2*NTC accumulated cosine terms !
47! !
48! CosW_avg Current time-average window COS(omega(k)*t). !
49! CosW_sum Time-accumulated COS(omega(k)*t). !
50! SinW_avg Current time-average window SIN(omega(k)*t). !
51! SinW_sum Time-accumulated SIN(omega(k)*t). !
52! CosWCosW Time-accumulated COS(omega(k)*t)*COS(omega(l)*t). !
53! SinWSinW Time-accumulated SIN(omega(k)*t)*SIN(omega(l)*t). !
54! SinWCosW Time-accumulated SIN(omega(k)*t)*COS(omega(l)*t). !
55! !
56! ubar_detided Time-averaged and detided 2D u-momentum. !
57! ubar_tide Time-accumulated 2D u-momentum tide harmonics. !
58! vbar_detided Time-averaged and detided 2D v-momentum. !
59! vbar_tide Time-accumulated 2D v-momentum tide harmonics. !
60! zeta_detided Time-averaged and detided free-surface. !
61! zeta_tide Time-accumulated free-surface tide harmonics. !
62# ifdef SOLVE3D
63! t_detided Time-averaged and detided tracers (T,S). !
64! t_tide Time-accumulated 3D tracers (T,S) tide harmonics. !
65! u_detided Time-averaged and detided 3D u-momentum. !
66! u_tide Time-accumulated 3D u-momentum tide harmonics. !
67! v_detided Time-averaged and detided 3D v-momentum. !
68! v_tide Time-accumulated 3D v-momentum tide harmonics. !
69# endif
70! !
71# endif
72!=======================================================================
73!
74 USE mod_kinds
75!
76 implicit none
77!
78 PUBLIC :: allocate_tides
79 PUBLIC :: deallocate_tides
80 PUBLIC :: initialize_tides
81!
82!-----------------------------------------------------------------------
83! Define T_TIDES structure.
84!-----------------------------------------------------------------------
85!
87
88 real(r8), pointer :: tperiod(:)
89# if defined AVERAGES && defined AVERAGES_DETIDE
90 real(r8), pointer :: cosomega(:)
91 real(r8), pointer :: sinomega(:)
92 real(r8), pointer :: cosw_avg(:)
93 real(r8), pointer :: cosw_sum(:)
94 real(r8), pointer :: sinw_avg(:)
95 real(r8), pointer :: sinw_sum(:)
96 real(r8), pointer :: coswcosw(:,:)
97 real(r8), pointer :: sinwsinw(:,:)
98 real(r8), pointer :: sinwcosw(:,:)
99# endif
100# if defined SSH_TIDES
101 real(r8), pointer :: ssh_tamp(:,:,:)
102 real(r8), pointer :: ssh_tphase(:,:,:)
103# endif
104# if defined UV_TIDES
105 real(r8), pointer :: uv_tangle(:,:,:)
106 real(r8), pointer :: uv_tmajor(:,:,:)
107 real(r8), pointer :: uv_tminor(:,:,:)
108 real(r8), pointer :: uv_tphase(:,:,:)
109# endif
110# if defined AVERAGES && defined AVERAGES_DETIDE
111 real(r8), pointer :: ubar_detided(:,:)
112 real(r8), pointer :: ubar_tide(:,:,:)
113
114 real(r8), pointer :: vbar_detided(:,:)
115 real(r8), pointer :: vbar_tide(:,:,:)
116
117 real(r8), pointer :: zeta_detided(:,:)
118 real(r8), pointer :: zeta_tide(:,:,:)
119# ifdef SOLVE3D
120 real(r8), pointer :: t_detided(:,:,:,:)
121 real(r8), pointer :: t_tide(:,:,:,:,:)
122
123 real(r8), pointer :: u_detided(:,:,:)
124 real(r8), pointer :: u_tide(:,:,:,:)
125
126 real(r8), pointer :: v_detided(:,:,:)
127 real(r8), pointer :: v_tide(:,:,:,:)
128# endif
129# endif
130
131 END TYPE t_tides
132!
133 TYPE (t_tides), allocatable :: tides(:)
134!
135 CONTAINS
136!
137 SUBROUTINE allocate_tides (ng, LBi, UBi, LBj, UBj)
138!
139!=======================================================================
140! !
141! This routine allocates all variables in the module for all nested !
142! grids. !
143! !
144!=======================================================================
145!
146 USE mod_param
147 USE mod_parallel
148 USE mod_iounits
149 USE mod_ncparam
150 USE mod_netcdf
151# if defined PIO_LIB && defined DISTRIBUTE
153# endif
154 USE mod_scalars
155 USE mod_stepping
156!
157 USE strings_mod, ONLY : founderror
158!
159! Inported variable declarations.
160!
161 integer, intent(in) :: ng, lbi, ubi, lbj, ubj
162!
163! Local variable declarations.
164!
165 logical :: foundit
166
167 integer :: nfiles, vid, i, ifile, mg, nvatt, nvdim
168
169 real(r8) :: size2d
170!
171 character (len=*), parameter :: myfile = &
172 & __FILE__//", allocate_tides"
173
174# if defined PIO_LIB && defined DISTRIBUTE
175!
176 TYPE (var_desc_t) :: my_piovar
177# endif
178!
179!-----------------------------------------------------------------------
180! Allocate module variables.
181!-----------------------------------------------------------------------
182!
183! Inquire about the maximum number of tidal components. Notice that
184! currently we only support nested applications where the tidal
185! forcing is applied to the main coarser grid (RefineScale(ng)=0)
186! and the other grids get the tidal forcing from the contact areas.
187!
188 IF (lprocesstides(ng)) THEN
189 mtc=0
190 foundit=.false.
191 SELECT CASE (tide(ng)%IOtype)
192 CASE (io_nf90)
193 CALL netcdf_inq_var (ng, inlm, tide(ng)%name, &
194 & myvarname = vname(1,idtper), &
195 & searchvar = foundit, &
196 & varid = vid, &
197 & nvardim = nvdim, &
198 & nvaratt = nvatt)
199
200# if defined PIO_LIB && defined DISTRIBUTE
201 CASE (io_pio)
202 CALL pio_netcdf_inq_var (ng, inlm, tide(ng)%name, &
203 & myvarname = vname(1,idtper), &
204 & searchvar = foundit, &
205 & piovar = my_piovar, &
206 & nvardim = nvdim, &
207 & nvaratt = nvatt)
208# endif
209 END SELECT
210 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
211!
212! Set maximum number of tidal components. Allocate and initialize
213! TIDE I/O structure. Notice that in nested applications, all the
214! nested grids need to have the same number of tidal component.
215!
216 IF (foundit) THEN
217 mtc=max(mtc,var_dsize(1)) ! first dimension
218 DO mg=1,ngrids
219 ntc(mg)=var_dsize(1)
220 END DO
221 END IF
222 END IF
223!
224! Allocate structure.
225!
226 IF (ng.eq.1) allocate ( tides(ngrids) )
227!
228! Set horizontal array size.
229!
230 size2d=real((ubi-lbi+1)*(ubj-lbj+1),r8)
231!
232! Allocate tidal forcing variables.
233!
234 allocate ( tides(ng) % Tperiod(mtc) )
235 dmem(ng)=dmem(ng)+real(mtc,r8)
236
237# if defined SSH_TIDES
238 allocate ( tides(ng) % SSH_Tamp(lbi:ubi,lbj:ubj,mtc) )
239 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
240
241 allocate ( tides(ng) % SSH_Tphase(lbi:ubi,lbj:ubj,mtc) )
242 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
243# endif
244
245# if defined UV_TIDES
246 allocate ( tides(ng) % UV_Tangle(lbi:ubi,lbj:ubj,mtc) )
247 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
248
249 allocate ( tides(ng) % UV_Tmajor(lbi:ubi,lbj:ubj,mtc) )
250 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
251
252 allocate ( tides(ng) % UV_Tminor(lbi:ubi,lbj:ubj,mtc) )
253 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
254
255 allocate ( tides(ng) % UV_Tphase(lbi:ubi,lbj:ubj,mtc) )
256 dmem(ng)=dmem(ng)+real(mtc,r8)*size2d
257# endif
258
259# if defined AVERAGES && defined AVERAGES_DETIDE
260!
261! Allocate variables used for the least-squares detiding of
262! time-averaged fields.
263!
264 allocate ( tides(ng) % CosOmega(mtc) )
265 dmem(ng)=dmem(ng)+real(mtc,r8)
266
267 allocate ( tides(ng) % SinOmega(mtc) )
268 dmem(ng)=dmem(ng)+real(mtc,r8)
269
270 allocate ( tides(ng) % CosW_avg(mtc) )
271 dmem(ng)=dmem(ng)+real(mtc,r8)
272
273 allocate ( tides(ng) % CosW_sum(mtc) )
274 dmem(ng)=dmem(ng)+real(mtc,r8)
275
276 allocate ( tides(ng) % SinW_avg(mtc) )
277 dmem(ng)=dmem(ng)+real(mtc,r8)
278
279 allocate ( tides(ng) % SinW_sum(mtc) )
280 dmem(ng)=dmem(ng)+real(mtc,r8)
281
282 allocate ( tides(ng) % CosWCosW(mtc,mtc) )
283 dmem(ng)=dmem(ng)+real(mtc*mtc,r8)
284
285 allocate ( tides(ng) % SinWSinW(mtc,mtc) )
286 dmem(ng)=dmem(ng)+real(mtc*mtc,r8)
287
288 allocate ( tides(ng) % SinWCosW(mtc,mtc) )
289 dmem(ng)=dmem(ng)+real(mtc*mtc,r8)
290
291 IF (aout(idfsud,ng)) THEN
292 allocate ( tides(ng) % zeta_detided(lbi:ubi,lbj:ubj) )
293 dmem(ng)=dmem(ng)+size2d
294
295 allocate ( tides(ng) % zeta_tide(lbi:ubi,lbj:ubj,0:2*mtc) )
296 dmem(ng)=dmem(ng)+real(2*mtc+1,r8)*size2d
297 END IF
298
299 IF (aout(idu2dd,ng)) THEN
300 allocate ( tides(ng) % ubar_detided(lbi:ubi,lbj:ubj) )
301
302 allocate ( tides(ng) % ubar_tide(lbi:ubi,lbj:ubj,0:2*mtc) )
303 dmem(ng)=dmem(ng)+real(2*mtc+1,r8)*size2d
304 END IF
305
306 IF (aout(idv2dd,ng)) THEN
307 allocate ( tides(ng) % vbar_detided(lbi:ubi,lbj:ubj) )
308
309 allocate ( tides(ng) % vbar_tide(lbi:ubi,lbj:ubj,0:2*mtc) )
310 dmem(ng)=dmem(ng)+real(2*mtc+1,r8)*size2d
311 END IF
312
313# ifdef SOLVE3D
314 IF (aout(idu3dd,ng)) THEN
315 allocate ( tides(ng) % u_detided(lbi:ubi,lbj:ubj,n(ng)) )
316 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
317
318 allocate ( tides(ng) % u_tide(lbi:ubi,lbj:ubj,n(ng),0:2*mtc) )
319 dmem(ng)=dmem(ng)+real(n(ng)*(2*mtc+1),r8)*size2d
320 END IF
321
322 IF (aout(idv3dd,ng)) THEN
323 allocate ( tides(ng) % v_detided(lbi:ubi,lbj:ubj,n(ng)) )
324 dmem(ng)=dmem(ng)+real(n(ng),r8)*size2d
325
326 allocate ( tides(ng) % v_tide(lbi:ubi,lbj:ubj,n(ng),0:2*mtc) )
327 dmem(ng)=dmem(ng)+real(n(ng)*(2*mtc+1),r8)*size2d
328 END IF
329
330 IF (any(aout(idtrcd(:),ng))) THEN
331 allocate ( tides(ng) % t_detided(lbi:ubi,lbj:ubj,n(ng),nat) )
332 dmem(ng)=dmem(ng)+real(n(ng)*nat,r8)*size2d
333
334 allocate ( tides(ng) % t_tide(lbi:ubi,lbj:ubj,n(ng), &
335 & 0:2*mtc,nat) )
336 dmem(ng)=dmem(ng)+real(n(ng)*(2*mtc+1)*nat,r8)*size2d
337 END IF
338# endif
339# endif
340!
341 RETURN
342 END SUBROUTINE allocate_tides
343!
344 SUBROUTINE deallocate_tides (ng)
345!
346!=======================================================================
347! !
348! This routine allocates all variables in the module for all nested !
349! grids. !
350! !
351!=======================================================================
352!
353 USE mod_param
354 USE mod_ncparam
355
356# ifdef SUBOBJECT_DEALLOCATION
357!
358 USE destroy_mod, ONLY : destroy
359# endif
360!
361! Inported variable declarations.
362!
363 integer, intent(in) :: ng
364!
365! Local variable declarations.
366!
367 character (len=*), parameter :: myfile = &
368 & __FILE__//", deallocate_tides"
369
370
371# ifdef SUBOBJECT_DEALLOCATION
372!
373!-----------------------------------------------------------------------
374! Deallocate each variable in the derived-type T_TIDES structure
375! separately.
376!-----------------------------------------------------------------------
377!
378 IF (.not.destroy(ng, tides(ng)%Tperiod, myfile, &
379 & __line__, 'TIDES(ng)%Tperiod')) RETURN
380
381# if defined SSH_TIDES
382 IF (.not.destroy(ng, tides(ng)%SSH_Tamp, myfile, &
383 & __line__, 'TIDES(ng)%SSH_Tamp')) RETURN
384
385 IF (.not.destroy(ng, tides(ng)%SSH_Tphase, myfile, &
386 & __line__, 'TIDES(ng)%SSH_Tphase')) RETURN
387# endif
388
389# if defined UV_TIDES
390 IF (.not.destroy(ng, tides(ng)%UV_Tangle, myfile, &
391 & __line__, 'TIDES(ng)%UV_Tangle')) RETURN
392
393 IF (.not.destroy(ng, tides(ng)%UV_Tmajor, myfile, &
394 & __line__, 'TIDES(ng)%UV_Tmajor')) RETURN
395
396 IF (.not.destroy(ng, tides(ng)%UV_Tminor, myfile, &
397 & __line__, 'TIDES(ng)%UV_Tminor')) RETURN
398
399 IF (.not.destroy(ng, tides(ng)%UV_Tphase, myfile, &
400 & __line__, 'TIDES(ng)%UV_Tphase')) RETURN
401# endif
402
403# if defined AVERAGES && defined AVERAGES_DETIDE
404!
405! Allocate variables used for the least-squares detiding of
406! time-averaged fields.
407!
408 IF (.not.destroy(ng, tides(ng)%CosOmega, myfile, &
409 & __line__, 'TIDES(ng)%CosOmega')) RETURN
410
411 IF (.not.destroy(ng, tides(ng)%SinOmega, myfile, &
412 & __line__, 'TIDES(ng)%SinOmega')) RETURN
413
414 IF (.not.destroy(ng, tides(ng)%CosW_avg, myfile, &
415 & __line__, 'TIDES(ng)%CosW_avg')) RETURN
416
417 IF (.not.destroy(ng, tides(ng)%CosW_sum, myfile, &
418 & __line__, 'TIDES(ng)%CosW_sum')) RETURN
419
420 IF (.not.destroy(ng, tides(ng)%SinW_avg, myfile, &
421 & __line__, 'TIDES(ng)%SinW_avg')) RETURN
422
423 IF (.not.destroy(ng, tides(ng)%SinW_sum, myfile, &
424 & __line__, 'TIDES(ng)%SinW_sum')) RETURN
425
426 IF (.not.destroy(ng, tides(ng)%CosWCosW, myfile, &
427 & __line__, 'TIDES(ng)%CosWCosW')) RETURN
428
429 IF (.not.destroy(ng, tides(ng)%SinWSinW, myfile, &
430 & __line__, 'TIDES(ng)%SinWSinW')) RETURN
431
432 IF (.not.destroy(ng, tides(ng)%SinWCosW, myfile, &
433 & __line__, 'TIDES(ng)%SinWCosW')) RETURN
434
435 IF (aout(idfsud,ng)) THEN
436 IF (.not.destroy(ng, tides(ng)%zeta_detided, myfile, &
437 & __line__, 'TIDES(ng)%zeta_detided')) RETURN
438
439 IF (.not.destroy(ng, tides(ng)%zeta_tide, myfile, &
440 & __line__, 'TIDES(ng)%zeta_tide')) RETURN
441 END IF
442
443 IF (aout(idu2dd,ng)) THEN
444 IF (.not.destroy(ng, tides(ng)%ubar_detided, myfile, &
445 & __line__, 'TIDES(ng)%ubar_detided')) RETURN
446
447 IF (.not.destroy(ng, tides(ng)%ubar_tide, myfile, &
448 & __line__, 'TIDES(ng)%ubar_tide')) RETURN
449 END IF
450
451 IF (aout(idv2dd,ng)) THEN
452 IF (.not.destroy(ng, tides(ng)%vbar_detided, myfile, &
453 & __line__, 'TIDES(ng)%vbar_detided')) RETURN
454
455 IF (.not.destroy(ng, tides(ng)%vbar_tide, myfile, &
456 & __line__, 'TIDES(ng)%vbar_tide')) RETURN
457 END IF
458
459# ifdef SOLVE3D
460 IF (aout(idu3dd,ng)) THEN
461 IF (.not.destroy(ng, tides(ng)%u_detided, myfile, &
462 & __line__, 'TIDES(ng)%u_detided')) RETURN
463
464 IF (.not.destroy(ng, tides(ng)%u_tide, myfile, &
465 & __line__, 'TIDES(ng)%u_tide')) RETURN
466 END IF
467
468 IF (aout(idv3dd,ng)) THEN
469 IF (.not.destroy(ng, tides(ng)%v_detided, myfile, &
470 & __line__, 'TIDES(ng)%v_detided')) RETURN
471
472 IF (.not.destroy(ng, tides(ng)%v_tide, myfile, &
473 & __line__, 'TIDES(ng)%v_tide')) RETURN
474 END IF
475
476 IF (any(aout(idtrcd(:),ng))) THEN
477 IF (.not.destroy(ng, tides(ng)%t_detided, myfile, &
478 & __line__, 'TIDES(ng)%t_detided')) RETURN
479
480 IF (.not.destroy(ng, tides(ng)%t_tide, myfile, &
481 & __line__, 'TIDES(ng)%t_tide')) RETURN
482 END IF
483# endif
484# endif
485# endif
486!
487!-----------------------------------------------------------------------
488! Deallocate derived-type TIDES structure.
489!-----------------------------------------------------------------------
490!
491 IF (ng.eq.ngrids) THEN
492 IF (allocated(tides)) deallocate ( tides )
493 END IF
494!
495 RETURN
496 END SUBROUTINE deallocate_tides
497!
498 SUBROUTINE initialize_tides (ng, tile)
499!
500!=======================================================================
501! !
502! This routine initialize all variables in the module using first !
503! touch distribution policy. In shared-memory configuration, this !
504! operation actually performs propagation of the "shared arrays" !
505! across the cluster, unless another policy is specified to !
506! override the default. !
507! !
508!=======================================================================
509!
510 USE mod_param
511 USE mod_ncparam
512!
513! Imported variable declarations.
514!
515 integer, intent(in) :: ng, tile
516!
517! Local variable declarations.
518!
519 integer :: imin, imax, jmin, jmax
520 integer :: i, itide, itrc, j, jtide, k
521
522 real(r8), parameter :: inival = 0.0_r8
523
524# include "set_bounds.h"
525!
526! Set array initialization range.
527!
528# ifdef DISTRIBUTE
529 imin=bounds(ng)%LBi(tile)
530 imax=bounds(ng)%UBi(tile)
531 jmin=bounds(ng)%LBj(tile)
532 jmax=bounds(ng)%UBj(tile)
533# else
534 IF (domain(ng)%Western_Edge(tile)) THEN
535 imin=bounds(ng)%LBi(tile)
536 ELSE
537 imin=istr
538 END IF
539 IF (domain(ng)%Eastern_Edge(tile)) THEN
540 imax=bounds(ng)%UBi(tile)
541 ELSE
542 imax=iend
543 END IF
544 IF (domain(ng)%Southern_Edge(tile)) THEN
545 jmin=bounds(ng)%LBj(tile)
546 ELSE
547 jmin=jstr
548 END IF
549 IF (domain(ng)%Northern_Edge(tile)) THEN
550 jmax=bounds(ng)%UBj(tile)
551 ELSE
552 jmax=jend
553 END IF
554# endif
555!
556!-----------------------------------------------------------------------
557! Initialize module variables.
558!-----------------------------------------------------------------------
559!
560! Initialize tidal forcing variables.
561!
562 IF (domain(ng)%SouthWest_Test(tile)) THEN
563 DO itide=1,mtc
564 tides(ng) % Tperiod(itide) = inival
565 END DO
566 END IF
567
568 DO itide=1,mtc
569# if defined SSH_TIDES
570 DO j=jmin,jmax
571 DO i=imin,imax
572 tides(ng) % SSH_Tamp(i,j,itide) = inival
573 tides(ng) % SSH_Tphase(i,j,itide) = inival
574 END DO
575 END DO
576# endif
577# if defined UV_TIDES
578 DO j=jmin,jmax
579 DO i=imin,imax
580 tides(ng) % UV_Tangle(i,j,itide) = inival
581 tides(ng) % UV_Tmajor(i,j,itide) = inival
582 tides(ng) % UV_Tminor(i,j,itide) = inival
583 tides(ng) % UV_Tphase(i,j,itide) = inival
584 END DO
585 END DO
586# endif
587 END DO
588
589# if defined AVERAGES && defined AVERAGES_DETIDE
590!
591! Initialize cariables used for the least-squares detiding of
592! time-averaged fields.
593!
594 IF (domain(ng)%SouthWest_Test(tile)) THEN
595 DO jtide=1,mtc
596 tides(ng) % CosOmega(jtide) = inival
597 tides(ng) % SinOmega(jtide) = inival
598 tides(ng) % CosW_avg(jtide) = inival
599 tides(ng) % CosW_sum(jtide) = inival
600 tides(ng) % SinW_avg(jtide) = inival
601 tides(ng) % SinW_sum(jtide) = inival
602 DO itide=1,mtc
603 tides(ng) % CosWCosW(itide,jtide) = inival
604 tides(ng) % SinWSinW(itide,jtide) = inival
605 tides(ng) % SinWCosW(itide,jtide) = inival
606 END DO
607 END DO
608 END IF
609
610 IF (aout(idfsud,ng)) THEN
611 DO j=jmin,jmax
612 DO i=imin,imax
613 tides(ng) % zeta_detided(i,j) = inival
614 END DO
615 END DO
616 DO itide=0,2*mtc
617 DO j=jmin,jmax
618 DO i=imin,imax
619 tides(ng) % zeta_tide(i,j,itide) = inival
620 END DO
621 END DO
622 END DO
623 END IF
624
625 IF (aout(idu2dd,ng)) THEN
626 DO j=jmin,jmax
627 DO i=imin,imax
628 tides(ng) % ubar_detided(i,j) = inival
629 END DO
630 END DO
631 DO itide=0,2*mtc
632 DO j=jmin,jmax
633 DO i=imin,imax
634 tides(ng) % ubar_tide(i,j,itide) = inival
635 END DO
636 END DO
637 END DO
638 END IF
639
640 IF (aout(idv2dd,ng)) THEN
641 DO j=jmin,jmax
642 DO i=imin,imax
643 tides(ng) % vbar_detided(i,j) = inival
644 END DO
645 END DO
646 DO itide=0,2*mtc
647 DO j=jmin,jmax
648 DO i=imin,imax
649 tides(ng) % vbar_tide(i,j,itide) = inival
650 END DO
651 END DO
652 END DO
653 END IF
654
655# ifdef SOLVE3D
656 IF (aout(idu3dd,ng)) THEN
657 DO k=1,n(ng)
658 DO j=jmin,jmax
659 DO i=imin,imax
660 tides(ng) % u_detided(i,j,k) = inival
661 END DO
662 END DO
663 END DO
664 DO itide=0,2*mtc
665 DO k=1,n(ng)
666 DO j=jmin,jmax
667 DO i=imin,imax
668 tides(ng) % u_tide(i,j,k,itide) = inival
669 END DO
670 END DO
671 END DO
672 END DO
673 END IF
674
675 IF (aout(idv3dd,ng)) THEN
676 DO k=1,n(ng)
677 DO j=jmin,jmax
678 DO i=imin,imax
679 tides(ng) % v_detided(i,j,k) = inival
680 END DO
681 END DO
682 END DO
683 DO itide=0,2*mtc
684 DO k=1,n(ng)
685 DO j=jmin,jmax
686 DO i=imin,imax
687 tides(ng) % v_tide(i,j,k,itide) = inival
688 END DO
689 END DO
690 END DO
691 END DO
692 END IF
693
694 IF (any(aout(idtrcd(:),ng))) THEN
695 DO itrc=1,nat
696 DO k=1,n(ng)
697 DO j=jmin,jmax
698 DO i=imin,imax
699 tides(ng) % t_detided(i,j,k,itrc) = inival
700 END DO
701 END DO
702 END DO
703 END DO
704 DO itrc=1,nat
705 DO itide=0,2*mtc
706 DO k=1,n(ng)
707 DO j=jmin,jmax
708 DO i=imin,imax
709 tides(ng) % t_tide(i,j,k,itide,itrc) = inival
710 END DO
711 END DO
712 END DO
713 END DO
714 END DO
715 END IF
716# endif
717# endif
718!
719 RETURN
720 END SUBROUTINE initialize_tides
721#endif
722 END MODULE mod_tides
type(t_io), dimension(:), allocatable tide
integer, parameter r8
Definition mod_kinds.F:28
integer, parameter io_nf90
Definition mod_ncparam.F:95
integer idu3dd
integer idv3dd
integer idtper
integer, parameter io_pio
Definition mod_ncparam.F:96
integer idfsud
integer, dimension(:), allocatable idtrcd
integer idv2dd
character(len=maxlen), dimension(6, 0:nv) vname
integer idu2dd
logical, dimension(:,:), allocatable aout
integer, dimension(nvard) var_dsize
Definition mod_netcdf.F:177
subroutine, public netcdf_inq_var(ng, model, ncname, ncid, myvarname, searchvar, varid, nvardim, nvaratt)
integer nat
Definition mod_param.F:499
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
integer mtc
Definition mod_param.F:564
real(r8), dimension(:), allocatable dmem
Definition mod_param.F:137
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer ngrids
Definition mod_param.F:113
subroutine, public pio_netcdf_inq_var(ng, model, ncname, piofile, myvarname, searchvar, piovar, nvardim, nvaratt)
logical, dimension(:), allocatable lprocesstides
integer exit_flag
integer noerror
integer, dimension(:), allocatable ntc
subroutine, public deallocate_tides(ng)
Definition mod_tides.F:345
subroutine, public initialize_tides(ng, tile)
Definition mod_tides.F:499
type(t_tides), dimension(:), allocatable tides
Definition mod_tides.F:133
subroutine, public allocate_tides(ng, lbi, ubi, lbj, ubj)
Definition mod_tides.F:138
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52